Flow123d  release_2.2.0-914-gf1a3a4f
Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Attributes | List of all members
TransportDG< Model > Class Template Reference

Transport with dispersion implemented using discontinuous Galerkin method. More...

#include <transport_dg.hh>

Inheritance diagram for TransportDG< Model >:
Inheritance graph
[legend]
Collaboration diagram for TransportDG< Model >:
Collaboration graph
[legend]

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)
 
IdxIntget_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::Recordget_input_type ()
 Declare input record type for the equation TransportDG. More...
 
static const Input::Type::Selectionget_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
FEObjectsfeo
 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< VectorSeqDoubleoutput_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...
 

Detailed Description

template<class Model>
class TransportDG< Model >

Transport with dispersion implemented using discontinuous Galerkin method.

TransportDG implements the discontinuous Galerkin method for the transport and diffusion of substances. The concentration $ c_i ~[kg/m^3]$ of the i-th substance is governed by the advection-diffusion equation

\[ \partial_t c_i + \mathbf v\cdot\nabla c_i - \mathrm{div}(D\nabla c_i) = F \mbox{ in }\Omega^d, \]

where $\mathbf v$ is the fluid velocity and $\Omega^d$ the $d$-dimensional domain, respectively. The hydrodynamic dispersivity tensor $\mathbf D ~[m^2/s]$ is given by:

\[ \mathbf D = D_m\mathbf I + |\mathbf v|\left(\alpha_T\mathbf I + (\alpha_L-\alpha_T)\frac{\mathbf v\otimes\mathbf v}{|\mathbf v|^2}\right). \]

The molecular dispersivity $D_m~[m^2/s]$, as well as the longitudal and transversal dispersivity $\alpha_L,~\alpha_T~[m]$ are input parameters of the model.

For lower dimensions $d=1,2$ the advection-diffusion equation is multiplied by the fracture cross-cut $\delta^d~[m^{3-d}]$.

The boundary $\partial\Omega^d$ is divided into three disjoint parts $\Gamma^d_D\cup\Gamma^d_N\cup\Gamma^d_F$. We prescribe the following boundary conditions:

\begin{eqnarray*} c_i^d &= c_{iD}^d &\mbox{ on }\Gamma^d_D \mbox{ (Dirichlet)},\\ \mathbf D^d\nabla c_i^d\cdot\mathbf n &= 0 &\mbox{ on }\Gamma^d_N \mbox{ (Neumann)}, \end{eqnarray*}

The transfer of mass through fractures is described by the transmission conditions on $\Gamma^d_F$:

\[ -\mathbf D^d\nabla c_i^d\cdot\mathbf n = \sigma(c_i^d-c_i^{d-1}) + \left\{\begin{array}{cl}0 &\mbox{ if }\mathbf v^d\cdot\mathbf n\ge 0\\\mathbf v^d\cdot\mathbf n(c_i^{d-1}-c_i^d) & \mbox{ if }\mathbf v^d\cdot\mathbf n<0\end{array}\right.,\qquad F^{d-1} = (\sigma + |\mathbf v^d\cdot\mathbf n|)(c_i^d-c_i^{d-1}). \]

Here $\mathbf n$ stands for the unit outward normal vector to $\partial\Omega^d$. The coefficient $\sigma$ determines the transfer of mass through fractures due to diffusion.

Definition at line 131 of file transport_dg.hh.

Member Enumeration Documentation

template<class Model >
enum TransportDG::DGVariant
Enumerator
non_symmetric 
incomplete 
symmetric 

Definition at line 151 of file transport_dg.hh.

Constructor & Destructor Documentation

template<class Model >
TransportDG< Model >::TransportDG ( Mesh init_mesh,
const Input::Record  in_rec 
)

Constructor.

Parameters
init_meshcomputational mesh
in_recinput record

Definition at line 205 of file transport_dg.cc.

template<class Model >
TransportDG< Model >::~TransportDG ( )

Destructor.

Definition at line 343 of file transport_dg.cc.

Member Function Documentation

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::assemble_fluxes_boundary ( )
private

Assembles the fluxes on the boundary.

Definition at line 1039 of file transport_dg.cc.

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::assemble_fluxes_element_element ( )
private

Assembles the fluxes between elements of the same dimension.

Definition at line 907 of file transport_dg.cc.

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::assemble_fluxes_element_side ( )
private

Assembles the fluxes between elements of different dimensions.

Definition at line 1140 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::assemble_mass_matrix ( )
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.

Here is the caller graph for this function:

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::assemble_mass_matrix ( )
private

Assembles the mass matrix for the given dimension.

template<class Model >
void TransportDG< Model >::assemble_stiffness_matrix ( )
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.

Here is the caller graph for this function:

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::assemble_volume_integrals ( )
private

Assembles the volume integrals into the stiffness matrix.

Definition at line 782 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::calculate_concentration_matrix ( )

Definition at line 595 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::calculate_cumulative_balance ( )

Definition at line 661 of file transport_dg.cc.

Here is the caller graph for this function:

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::calculate_velocity ( const ElementFullIter cell,
std::vector< arma::vec3 > &  velocity,
FEValuesBase< dim, 3 > &  fv 
)
private

Calculates the velocity field on a given dim dimensional cell.

Parameters
cellThe cell.
velocityThe computed velocity field (at quadrature points).
fvThe FEValues class providing the quadrature points and the shape functions for velocity.

Here is the caller graph for this function:

template<class Model >
template<unsigned int dim>
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.

template<class Model >
Model::ModelEqData& TransportDG< Model >::data ( )
inlineprivate

Definition at line 227 of file transport_dg.hh.

template<class Model >
bool TransportDG< Model >::evaluate_time_constraint ( double &  time_constraint)
inline

Definition at line 184 of file transport_dg.hh.

template<class Model >
double** TransportDG< Model >::get_concentration_matrix ( )
inline

Definition at line 209 of file transport_dg.hh.

template<class Model >
const Selection & TransportDG< Model >::get_dg_variant_selection_input_type ( )
static

Input type for the DG variant selection.

Definition at line 48 of file transport_dg.cc.

template<class Model >
const Record & TransportDG< Model >::get_input_type ( )
static

Declare input record type for the equation TransportDG.

Definition at line 74 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::get_par_info ( IdxInt *&  el_4_loc,
Distribution *&  el_ds 
)

Definition at line 1683 of file transport_dg.cc.

template<class Model >
IdxInt * TransportDG< Model >::get_row_4_el ( )

Definition at line 1734 of file transport_dg.cc.

template<class Model >
const Vec& TransportDG< Model >::get_solution ( unsigned int  sbi)
inline

Definition at line 206 of file transport_dg.hh.

Here is the caller graph for this function:

template<class Model >
void TransportDG< Model >::initialize ( )
override

Definition at line 239 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::output_data ( void  )

Postprocesses the solution and writes to output file.

Definition at line 635 of file transport_dg.cc.

Here is the caller graph for this function:

template<class Model >
void TransportDG< Model >::output_vector_gather ( void  )
private

Definition at line 383 of file transport_dg.cc.

Here is the caller graph for this function:

template<class Model >
void TransportDG< Model >::preallocate ( )
private

Definition at line 433 of file transport_dg.cc.

Here is the caller graph for this function:

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::prepare_initial_condition ( )
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.

template<class Model >
void TransportDG< Model >::set_boundary_conditions ( )
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.

Here is the caller graph for this function:

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::set_boundary_conditions ( )
private

Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions for a given space dimension.

template<class Model >
void TransportDG< Model >::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 
)
private

Sets up parameters of the DG method on a given boundary edge.

Assumption is that the edge consists of only 1 side.

Parameters
sideThe boundary side.
K_sizeSize of vector of tensors K.
KDispersivity tensor.
ad_vectorAdvection vector.
normal_vectorNormal vector (assumed constant along the edge).
alphaPenalty parameter that influences the continuity of the solution (large value=more continuity).
gammaComputed penalty parameters.

Definition at line 1572 of file transport_dg.cc.

Here is the caller graph for this function:

template<class Model >
void TransportDG< Model >::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 
)
private

Calculates the dispersivity (diffusivity) tensor from the velocity field.

Parameters
KThe computed dispersivity tensor.
velocityThe velocity field (at quadrature points).
DmMolecular diffusivities.
alphaLLongitudal dispersivities.
alphaTTransversal dispersivities.
porosityPorosities.
cross_cutCross-cuts of higher dimension. Sets up some parameters of the DG method for two sides of an edge.
edgThe edge.
s1Side 1.
s2Side 2.
K_sizeSize of vector of tensors K.
K1Dispersivity tensors on side s1 (in quadrature points).
K2Dispersivity tensors on side s2 (in quadrature points).
normal_vectorNormal vector to side 0 of the neighbour (assumed constant along the side).
alpha1,alpha2Penalty parameter that influences the continuity of the solution (large value=more continuity).
gammaComputed penalty parameters.
omegaComputed weights.
transport_fluxComputed flux from side s1 to side s2.

Definition at line 1465 of file transport_dg.cc.

Here is the caller graph for this function:

template<class Model >
void TransportDG< Model >::set_initial_condition ( )
private

Sets the initial condition.

Definition at line 1607 of file transport_dg.cc.

Here is the caller graph for this function:

template<class Model >
void TransportDG< Model >::set_sources ( )
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.

Here is the caller graph for this function:

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::set_sources ( )
private

Assembles the right hand side vector due to volume sources.

template<class Model >
void TransportDG< Model >::update_after_reactions ( bool  solution_changed)

Definition at line 1691 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::update_solution ( void  )
override

Computes the solution in one time instant.

Definition at line 464 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::zero_time_step ( )
override

Initialize solution in the zero time.

Definition at line 400 of file transport_dg.cc.

Member Data Documentation

template<class Model >
vector<vector<arma::vec3> > TransportDG< Model >::ad_coef
private

Advection coefficients.

Definition at line 485 of file transport_dg.hh.

template<class Model >
vector<vector<vector<arma::vec3> > > TransportDG< Model >::ad_coef_edg
private

Advection coefficients on edges.

Definition at line 489 of file transport_dg.hh.

template<class Model >
bool TransportDG< Model >::allocation_done
private

Indicates whether matrices have been preallocated.

Definition at line 502 of file transport_dg.hh.

template<class Model >
EqData TransportDG< Model >::data_
private

Field data for model parameters.

Definition at line 405 of file transport_dg.hh.

template<class Model >
unsigned int TransportDG< Model >::dg_order
private

Polynomial order of finite elements.

Definition at line 423 of file transport_dg.hh.

template<class Model >
int TransportDG< Model >::dg_variant
private

DG variant ((non-)symmetric/incomplete.

Definition at line 420 of file transport_dg.hh.

template<class Model >
vector<vector<arma::mat33> > TransportDG< Model >::dif_coef
private

Diffusion coefficients.

Definition at line 487 of file transport_dg.hh.

template<class Model >
vector<vector<vector<arma::mat33> > > TransportDG< Model >::dif_coef_edg
private

Diffusion coefficients on edges.

Definition at line 491 of file transport_dg.hh.

template<class Model >
FEObjects* TransportDG< Model >::feo
private

Finite element objects.

Definition at line 414 of file transport_dg.hh.

template<class Model >
std::vector<std::vector<double> > TransportDG< Model >::gamma
private

Penalty parameters.

Definition at line 417 of file transport_dg.hh.

template<class Model >
Input::Record TransportDG< Model >::input_rec
private

Record with input specification.

Definition at line 469 of file transport_dg.hh.

template<class Model >
LinSys** TransportDG< Model >::ls
private

Linear algebra system for the transport equation.

Definition at line 448 of file transport_dg.hh.

template<class Model >
LinSys** TransportDG< Model >::ls_dt
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.

template<class Model >
std::vector<Mat> TransportDG< Model >::mass_matrix
private

The mass matrix.

Definition at line 439 of file transport_dg.hh.

template<class Model >
std::vector<Vec> TransportDG< Model >::mass_vec
private

Mass from previous time instant (necessary when coefficients of mass matrix change in time).

Definition at line 442 of file transport_dg.hh.

template<class Model >
vector<double> TransportDG< Model >::mm_coef
private

Mass matrix coefficients.

Definition at line 479 of file transport_dg.hh.

template<class Model >
vector<VectorSeqDouble> TransportDG< Model >::output_vec
private

Array for storing the output solution data.

Vector of solution data.

Definition at line 466 of file transport_dg.hh.

template<class Model >
const int TransportDG< Model >::registrar
staticprivate
Initial value:
=
Input::register_class< TransportDG<Model>, Mesh &, const Input::Record>(std::string(Model::ModelEqData::name()) + "_DG") +

Registrar of class to factory.

Definition at line 225 of file transport_dg.hh.

template<class Model >
vector<vector<double> > TransportDG< Model >::ret_coef
private

Retardation coefficient due to sorption.

Definition at line 481 of file transport_dg.hh.

template<class Model >
vector<double> TransportDG< Model >::ret_sources
private

Temporary values of increments due to retardation (e.g. sorption)

Definition at line 483 of file transport_dg.hh.

template<class Model >
vector<double> TransportDG< Model >::ret_sources_prev
private

Definition at line 483 of file transport_dg.hh.

template<class Model >
std::vector<Vec> TransportDG< Model >::ret_vec
private

Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).

Definition at line 445 of file transport_dg.hh.

template<class Model >
std::vector<Vec> TransportDG< Model >::rhs
private

Vector of right hand side.

Definition at line 433 of file transport_dg.hh.

template<class Model >
double** TransportDG< Model >::solution_elem_
private

Element averages of solution (the array is passed to reactions in operator splitting).

Definition at line 454 of file transport_dg.hh.

template<class Model >
std::vector<Mat> TransportDG< Model >::stiffness_matrix
private

The stiffness matrix.

Definition at line 436 of file transport_dg.hh.


The documentation for this class was generated from the following files: