Flow123d
release_2.2.0-914-gf1a3a4f
|
Transport with dispersion implemented using discontinuous Galerkin method. More...
#include <transport_dg.hh>
Classes | |
class | EqData |
Public Types | |
enum | DGVariant { non_symmetric = -1, incomplete = 0, symmetric = 1 } |
Public Member Functions | |
TransportDG (Mesh &init_mesh, const Input::Record in_rec) | |
Constructor. More... | |
void | zero_time_step () override |
Initialize solution in the zero time. More... | |
bool | evaluate_time_constraint (double &time_constraint) |
void | update_solution () override |
Computes the solution in one time instant. More... | |
void | output_data () |
Postprocesses the solution and writes to output file. More... | |
~TransportDG () | |
Destructor. More... | |
void | initialize () override |
void | calculate_cumulative_balance () |
const Vec & | get_solution (unsigned int sbi) |
double ** | get_concentration_matrix () |
void | calculate_concentration_matrix () |
void | update_after_reactions (bool solution_changed) |
void | get_par_info (IdxInt *&el_4_loc, Distribution *&el_ds) |
IdxInt * | get_row_4_el () |
template<unsigned int dim> | |
void | calculate_velocity (const typename DOFHandlerBase::CellIterator &cell, vector< arma::vec3 > &velocity, FEValuesBase< dim, 3 > &fv) |
Static Public Member Functions | |
static const Input::Type::Record & | get_input_type () |
Declare input record type for the equation TransportDG. More... | |
static const Input::Type::Selection & | get_dg_variant_selection_input_type () |
Input type for the DG variant selection. More... | |
Private Member Functions | |
Model::ModelEqData & | data () |
void | output_vector_gather () |
void | preallocate () |
void | assemble_mass_matrix () |
Assembles the mass matrix. More... | |
template<unsigned int dim> | |
void | assemble_mass_matrix () |
Assembles the mass matrix for the given dimension. More... | |
void | assemble_stiffness_matrix () |
Assembles the stiffness matrix. More... | |
template<unsigned int dim> | |
void | assemble_volume_integrals () |
Assembles the volume integrals into the stiffness matrix. More... | |
void | set_sources () |
Assembles the right hand side due to volume sources. More... | |
template<unsigned int dim> | |
void | set_sources () |
Assembles the right hand side vector due to volume sources. More... | |
template<unsigned int dim> | |
void | assemble_fluxes_boundary () |
Assembles the fluxes on the boundary. More... | |
template<unsigned int dim> | |
void | assemble_fluxes_element_element () |
Assembles the fluxes between elements of the same dimension. More... | |
template<unsigned int dim> | |
void | assemble_fluxes_element_side () |
Assembles the fluxes between elements of different dimensions. More... | |
void | set_boundary_conditions () |
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions. More... | |
template<unsigned int dim> | |
void | set_boundary_conditions () |
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions for a given space dimension. More... | |
template<unsigned int dim> | |
void | calculate_velocity (const ElementFullIter &cell, std::vector< arma::vec3 > &velocity, FEValuesBase< dim, 3 > &fv) |
Calculates the velocity field on a given dim dimensional cell. More... | |
void | set_DG_parameters_edge (const Edge &edg, const int s1, const int s2, const int K_size, const std::vector< arma::mat33 > &K1, const std::vector< arma::mat33 > &K2, const std::vector< double > &fluxes, const arma::vec3 &normal_vector, const double alpha1, const double alpha2, double &gamma, double *omega, double &transport_flux) |
Calculates the dispersivity (diffusivity) tensor from the velocity field. More... | |
void | set_DG_parameters_boundary (const SideIter side, const int K_size, const std::vector< arma::mat33 > &K, const double flux, const arma::vec3 &normal_vector, const double alpha, double &gamma) |
Sets up parameters of the DG method on a given boundary edge. More... | |
void | set_initial_condition () |
Sets the initial condition. More... | |
template<unsigned int dim> | |
void | prepare_initial_condition () |
Assembles the auxiliary linear system to calculate the initial solution as L^2-projection of the prescribed initial condition. More... | |
Private Attributes | |
Physical parameters | |
EqData | data_ |
Field data for model parameters. More... | |
Parameters of the numerical method | |
FEObjects * | feo |
Finite element objects. More... | |
std::vector< std::vector< double > > | gamma |
Penalty parameters. More... | |
int | dg_variant |
DG variant ((non-)symmetric/incomplete. More... | |
unsigned int | dg_order |
Polynomial order of finite elements. More... | |
Solution of algebraic system | |
std::vector< Vec > | rhs |
Vector of right hand side. More... | |
std::vector< Mat > | stiffness_matrix |
The stiffness matrix. More... | |
std::vector< Mat > | mass_matrix |
The mass matrix. More... | |
std::vector< Vec > | mass_vec |
Mass from previous time instant (necessary when coefficients of mass matrix change in time). More... | |
std::vector< Vec > | ret_vec |
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption). More... | |
LinSys ** | ls |
Linear algebra system for the transport equation. More... | |
LinSys ** | ls_dt |
Linear algebra system for the time derivative (actually it is used only for handling the matrix structures). More... | |
double ** | solution_elem_ |
Element averages of solution (the array is passed to reactions in operator splitting). More... | |
Output to file | |
vector< VectorSeqDouble > | output_vec |
Array for storing the output solution data. More... | |
Input::Record | input_rec |
Record with input specification. More... | |
Auxiliary fields used during assembly | |
vector< double > | mm_coef |
Mass matrix coefficients. More... | |
vector< vector< double > > | ret_coef |
Retardation coefficient due to sorption. More... | |
vector< double > | ret_sources |
Temporary values of increments due to retardation (e.g. sorption) More... | |
vector< double > | ret_sources_prev |
vector< vector< arma::vec3 > > | ad_coef |
Advection coefficients. More... | |
vector< vector< arma::mat33 > > | dif_coef |
Diffusion coefficients. More... | |
vector< vector< vector< arma::vec3 > > > | ad_coef_edg |
Advection coefficients on edges. More... | |
vector< vector< vector< arma::mat33 > > > | dif_coef_edg |
Diffusion coefficients on edges. More... | |
Other | |
bool | allocation_done |
Indicates whether matrices have been preallocated. More... | |
Static Private Attributes | |
static const int | registrar |
Registrar of class to factory. More... | |
Transport with dispersion implemented using discontinuous Galerkin method.
TransportDG implements the discontinuous Galerkin method for the transport and diffusion of substances. The concentration of the i-th substance is governed by the advection-diffusion equation
where is the fluid velocity and the -dimensional domain, respectively. The hydrodynamic dispersivity tensor is given by:
The molecular dispersivity , as well as the longitudal and transversal dispersivity are input parameters of the model.
For lower dimensions the advection-diffusion equation is multiplied by the fracture cross-cut .
The boundary is divided into three disjoint parts . We prescribe the following boundary conditions:
The transfer of mass through fractures is described by the transmission conditions on :
Here stands for the unit outward normal vector to . The coefficient determines the transfer of mass through fractures due to diffusion.
Definition at line 131 of file transport_dg.hh.
enum TransportDG::DGVariant |
Enumerator | |
---|---|
non_symmetric | |
incomplete | |
symmetric |
Definition at line 151 of file transport_dg.hh.
TransportDG< Model >::TransportDG | ( | Mesh & | init_mesh, |
const Input::Record | in_rec | ||
) |
Constructor.
init_mesh | computational mesh |
in_rec | input record |
Definition at line 205 of file transport_dg.cc.
TransportDG< Model >::~TransportDG | ( | ) |
Destructor.
Definition at line 343 of file transport_dg.cc.
|
private |
Assembles the fluxes on the boundary.
Definition at line 1039 of file transport_dg.cc.
|
private |
Assembles the fluxes between elements of the same dimension.
Definition at line 907 of file transport_dg.cc.
|
private |
Assembles the fluxes between elements of different dimensions.
Definition at line 1140 of file transport_dg.cc.
|
private |
Assembles the mass matrix.
The routine just calls templated method assemble_mass_matrix() for each space dimension.
Definition at line 680 of file transport_dg.cc.
|
private |
Assembles the mass matrix for the given dimension.
|
private |
Assembles the stiffness matrix.
This routine just calls assemble_volume_integrals(), assemble_fluxes_boundary(), assemble_fluxes_element_element() and assemble_fluxes_element_side() for each space dimension.
Definition at line 749 of file transport_dg.cc.
|
private |
Assembles the volume integrals into the stiffness matrix.
Definition at line 782 of file transport_dg.cc.
void TransportDG< Model >::calculate_concentration_matrix | ( | ) |
Definition at line 595 of file transport_dg.cc.
void TransportDG< Model >::calculate_cumulative_balance | ( | ) |
|
private |
Calculates the velocity field on a given dim
dimensional cell.
cell | The cell. |
velocity | The computed velocity field (at quadrature points). |
fv | The FEValues class providing the quadrature points and the shape functions for velocity. |
void TransportDG< Model >::calculate_velocity | ( | const typename DOFHandlerBase::CellIterator & | cell, |
vector< arma::vec3 > & | velocity, | ||
FEValuesBase< dim, 3 > & | fv | ||
) |
Definition at line 1431 of file transport_dg.cc.
|
inlineprivate |
Definition at line 227 of file transport_dg.hh.
|
inline |
Definition at line 184 of file transport_dg.hh.
|
inline |
Definition at line 209 of file transport_dg.hh.
|
static |
Input type for the DG variant selection.
Definition at line 48 of file transport_dg.cc.
|
static |
Declare input record type for the equation TransportDG.
Definition at line 74 of file transport_dg.cc.
void TransportDG< Model >::get_par_info | ( | IdxInt *& | el_4_loc, |
Distribution *& | el_ds | ||
) |
Definition at line 1683 of file transport_dg.cc.
IdxInt * TransportDG< Model >::get_row_4_el | ( | ) |
Definition at line 1734 of file transport_dg.cc.
|
inline |
|
override |
Definition at line 239 of file transport_dg.cc.
void TransportDG< Model >::output_data | ( | void | ) |
Postprocesses the solution and writes to output file.
Definition at line 635 of file transport_dg.cc.
|
private |
|
private |
|
private |
Assembles the auxiliary linear system to calculate the initial solution as L^2-projection of the prescribed initial condition.
Definition at line 1632 of file transport_dg.cc.
|
private |
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
The routine just calls templated method set_boundary_condition() for each space dimension.
Definition at line 1259 of file transport_dg.cc.
|
private |
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions for a given space dimension.
|
private |
Sets up parameters of the DG method on a given boundary edge.
Assumption is that the edge consists of only 1 side.
side | The boundary side. |
K_size | Size of vector of tensors K. |
K | Dispersivity tensor. |
ad_vector | Advection vector. |
normal_vector | Normal vector (assumed constant along the edge). |
alpha | Penalty parameter that influences the continuity of the solution (large value=more continuity). |
gamma | Computed penalty parameters. |
Definition at line 1572 of file transport_dg.cc.
|
private |
Calculates the dispersivity (diffusivity) tensor from the velocity field.
K | The computed dispersivity tensor. |
velocity | The velocity field (at quadrature points). |
Dm | Molecular diffusivities. |
alphaL | Longitudal dispersivities. |
alphaT | Transversal dispersivities. |
porosity | Porosities. |
cross_cut | Cross-cuts of higher dimension. Sets up some parameters of the DG method for two sides of an edge. |
edg | The edge. |
s1 | Side 1. |
s2 | Side 2. |
K_size | Size of vector of tensors K. |
K1 | Dispersivity tensors on side s1 (in quadrature points). |
K2 | Dispersivity tensors on side s2 (in quadrature points). |
normal_vector | Normal vector to side 0 of the neighbour (assumed constant along the side). |
alpha1,alpha2 | Penalty parameter that influences the continuity of the solution (large value=more continuity). |
gamma | Computed penalty parameters. |
omega | Computed weights. |
transport_flux | Computed flux from side s1 to side s2. |
Definition at line 1465 of file transport_dg.cc.
|
private |
Sets the initial condition.
Definition at line 1607 of file transport_dg.cc.
|
private |
Assembles the right hand side due to volume sources.
This method just calls set_sources() for each space dimension.
Definition at line 836 of file transport_dg.cc.
|
private |
Assembles the right hand side vector due to volume sources.
void TransportDG< Model >::update_after_reactions | ( | bool | solution_changed | ) |
Definition at line 1691 of file transport_dg.cc.
|
override |
Computes the solution in one time instant.
Definition at line 464 of file transport_dg.cc.
|
override |
Initialize solution in the zero time.
Definition at line 400 of file transport_dg.cc.
|
private |
Advection coefficients.
Definition at line 485 of file transport_dg.hh.
|
private |
Advection coefficients on edges.
Definition at line 489 of file transport_dg.hh.
|
private |
Indicates whether matrices have been preallocated.
Definition at line 502 of file transport_dg.hh.
|
private |
Field data for model parameters.
Definition at line 405 of file transport_dg.hh.
|
private |
Polynomial order of finite elements.
Definition at line 423 of file transport_dg.hh.
|
private |
DG variant ((non-)symmetric/incomplete.
Definition at line 420 of file transport_dg.hh.
|
private |
Diffusion coefficients.
Definition at line 487 of file transport_dg.hh.
|
private |
Diffusion coefficients on edges.
Definition at line 491 of file transport_dg.hh.
|
private |
Finite element objects.
Definition at line 414 of file transport_dg.hh.
|
private |
Penalty parameters.
Definition at line 417 of file transport_dg.hh.
|
private |
Record with input specification.
Definition at line 469 of file transport_dg.hh.
|
private |
Linear algebra system for the transport equation.
Definition at line 448 of file transport_dg.hh.
|
private |
Linear algebra system for the time derivative (actually it is used only for handling the matrix structures).
Definition at line 451 of file transport_dg.hh.
|
private |
The mass matrix.
Definition at line 439 of file transport_dg.hh.
|
private |
Mass from previous time instant (necessary when coefficients of mass matrix change in time).
Definition at line 442 of file transport_dg.hh.
|
private |
Mass matrix coefficients.
Definition at line 479 of file transport_dg.hh.
|
private |
Array for storing the output solution data.
Vector of solution data.
Definition at line 466 of file transport_dg.hh.
|
staticprivate |
Registrar of class to factory.
Definition at line 225 of file transport_dg.hh.
|
private |
Retardation coefficient due to sorption.
Definition at line 481 of file transport_dg.hh.
|
private |
Temporary values of increments due to retardation (e.g. sorption)
Definition at line 483 of file transport_dg.hh.
|
private |
Definition at line 483 of file transport_dg.hh.
|
private |
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
Definition at line 445 of file transport_dg.hh.
|
private |
Vector of right hand side.
Definition at line 433 of file transport_dg.hh.
|
private |
Element averages of solution (the array is passed to reactions in operator splitting).
Definition at line 454 of file transport_dg.hh.
|
private |
The stiffness matrix.
Definition at line 436 of file transport_dg.hh.