Flow123d  JS_before_hm-1001-gfa0c761
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 }
 
template<unsigned int dim>
using MassAssemblyDim = MassAssemblyDG< dim, Model >
 
template<unsigned int dim>
using StiffnessAssemblyDim = StiffnessAssemblyDG< dim, Model >
 
template<unsigned int dim>
using SourcesAssemblyDim = SourcesAssemblyDG< dim, Model >
 
template<unsigned int dim>
using BdrConditionAssemblyDim = BdrConditionAssemblyDG< dim, Model >
 
template<unsigned int dim>
using InitConditionAssemblyDim = InitConditionAssemblyDG< dim, Model >
 
typedef std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > FieldFEScalarVec
 
- Public Types inherited from Model< spacedim, Value >
typedef FieldAlgorithmBase< spacedim, ValueFieldBaseType
 
typedef std::shared_ptr< FieldBaseTypeFieldBasePtr
 

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 &)
 
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 () override
 Destructor. More...
 
void initialize () override
 
void calculate_cumulative_balance ()
 
Vec get_component_vec (unsigned int sbi)
 Return PETSc vector with solution for sbi-th component. More...
 
FieldFEScalarVecget_p0_interpolation ()
 Getter for P0 interpolation by FieldFE. More...
 
void compute_p0_interpolation ()
 Compute P0 interpolation of the solution (used in reaction term). More...
 
void update_after_reactions (bool solution_changed)
 
void get_par_info (LongIdx *&el_4_loc, Distribution *&el_ds)
 
LongIdxget_row_4_el ()
 
std::shared_ptr< Balancebalance () const
 Access to balance object of Model. More...
 
const vector< unsigned int > subst_idx () const
 Return vector of substances indices. More...
 
Model::ModelEqData & data ()
 

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...
 
- Static Public Member Functions inherited from Model< spacedim, Value >
template<typename Fn , class... InputFields>
static auto create (Fn *fn, InputFields &&...inputs) -> decltype(auto)
 
template<typename Function , typename Tuple , size_t... I>
static auto call_create (Function f, Tuple t, std::index_sequence< I... >)
 
template<typename Fn , class... InputFields>
static auto create_multi (Fn *fn, InputFields &&...inputs) -> decltype(auto)
 

Private Member Functions

void preallocate ()
 
void set_initial_condition ()
 Calculates the dispersivity (diffusivity) tensor from the velocity field. More...
 
void output_region_statistics ()
 

Private Attributes

Physical parameters
std::shared_ptr< EqDatadata_
 Field data for model parameters. 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...
 
Output to file
vector< VectorMPIoutput_vec
 Array for storing the output solution data. More...
 
Input::Record input_rec
 Record with input specification. More...
 
ofstream reg_stat_stream
 
Auxiliary fields used during assembly
vector< double > ret_sources
 Temporary values of increments due to retardation (e.g. sorption) More...
 
vector< double > ret_sources_prev
 
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 132 of file transport_dg.hh.

Member Typedef Documentation

template<class Model>
template<unsigned int dim>
using TransportDG< Model >::BdrConditionAssemblyDim = BdrConditionAssemblyDG<dim, Model>

Definition at line 139 of file transport_dg.hh.

template<class Model>
typedef std::vector<std::shared_ptr<FieldFE< 3, FieldValue<3>::Scalar> > > TransportDG< Model >::FieldFEScalarVec

Definition at line 142 of file transport_dg.hh.

template<class Model>
template<unsigned int dim>
using TransportDG< Model >::InitConditionAssemblyDim = InitConditionAssemblyDG<dim, Model>

Definition at line 140 of file transport_dg.hh.

template<class Model>
template<unsigned int dim>
using TransportDG< Model >::MassAssemblyDim = MassAssemblyDG<dim, Model>

Definition at line 136 of file transport_dg.hh.

template<class Model>
template<unsigned int dim>
using TransportDG< Model >::SourcesAssemblyDim = SourcesAssemblyDG<dim, Model>

Definition at line 138 of file transport_dg.hh.

template<class Model>
template<unsigned int dim>
using TransportDG< Model >::StiffnessAssemblyDim = StiffnessAssemblyDG<dim, Model>

Definition at line 137 of file transport_dg.hh.

Member Enumeration Documentation

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

Definition at line 236 of file transport_dg.hh.

Constructor & Destructor Documentation

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

Constructor.

Parameters
init_meshcomputational mesh
in_recinput record

Definition at line 198 of file transport_dg.cc.

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

Destructor.

Definition at line 356 of file transport_dg.cc.

Member Function Documentation

template<class Model>
std::shared_ptr<Balance> TransportDG< Model >::balance ( ) const
inline

Access to balance object of Model.

Definition at line 309 of file transport_dg.hh.

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

Definition at line 659 of file transport_dg.cc.

Here is the caller graph for this function:

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

Compute P0 interpolation of the solution (used in reaction term).

Definition at line 607 of file transport_dg.cc.

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

Definition at line 319 of file transport_dg.hh.

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

Definition at line 269 of file transport_dg.hh.

template<class Model>
Vec TransportDG< Model >::get_component_vec ( unsigned int  sbi)
inline

Return PETSc vector with solution for sbi-th component.

Definition at line 292 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 52 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 78 of file transport_dg.cc.

template<class Model>
FieldFEScalarVec& TransportDG< Model >::get_p0_interpolation ( )
inline

Getter for P0 interpolation by FieldFE.

Definition at line 296 of file transport_dg.hh.

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

Definition at line 701 of file transport_dg.cc.

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

Definition at line 740 of file transport_dg.cc.

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

Definition at line 243 of file transport_dg.cc.

Here is the caller graph for this function:

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

Postprocesses the solution and writes to output file.

Definition at line 634 of file transport_dg.cc.

Here is the caller graph for this function:

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

Definition at line 431 of file transport_dg.cc.

Here is the caller graph for this function:

template<class Model >
void TransportDG< Model >::set_initial_condition ( )
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 the initial condition.

This method just calls AssemblyDG::prepare_initial_condition() for each elements.

Definition at line 680 of file transport_dg.cc.

Here is the caller graph for this function:

template<class Model>
const vector<unsigned int> TransportDG< Model >::subst_idx ( ) const
inline

Return vector of substances indices.

Definition at line 314 of file transport_dg.hh.

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

Definition at line 709 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 470 of file transport_dg.cc.

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

Initialize solution in the zero time.

Definition at line 398 of file transport_dg.cc.

Member Data Documentation

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

Indicates whether matrices have been preallocated.

Definition at line 415 of file transport_dg.hh.

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

Field data for model parameters.

Definition at line 359 of file transport_dg.hh.

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

Record with input specification.

Definition at line 392 of file transport_dg.hh.

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

The mass matrix.

Definition at line 375 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 378 of file transport_dg.hh.

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

Array for storing the output solution data.

Vector of solution data.

Definition at line 389 of file transport_dg.hh.

template<class Model>
ofstream TransportDG< Model >::reg_stat_stream
private

Definition at line 394 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 323 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 404 of file transport_dg.hh.

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

Definition at line 404 of file transport_dg.hh.

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

Vector of right hand side.

Definition at line 369 of file transport_dg.hh.

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

The stiffness matrix.

Definition at line 372 of file transport_dg.hh.


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