Flow123d
3.9.0-72cb42cb8
|
Transport with dispersion implemented using discontinuous Galerkin method. More...
#include <transport_dg.hh>
Classes | |
class | EqData |
class | EqFields |
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 > |
template<unsigned int dim> | |
using | InitProjectionAssemblyDim = InitProjectionAssemblyDG< dim, Model > |
typedef std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > | FieldFEScalarVec |
Public Types inherited from Model< spacedim, Value > | |
typedef FieldAlgorithmBase< spacedim, Value > | FieldBaseType |
typedef std::shared_ptr< FieldBaseType > | FieldBasePtr |
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... | |
FieldFEScalarVec & | get_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) |
std::shared_ptr< Balance > | balance () const |
Access to balance object of Model. More... | |
Model::ModelEqFields & | eq_fields () |
Model::ModelEqData & | eq_data () |
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... | |
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< EqFields > | eq_fields_ |
Fields for model parameters. More... | |
std::shared_ptr< EqData > | eq_data_ |
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 | |
Input::Record | input_rec |
Array for storing the output solution data. 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... | |
bool | init_projection |
GenericAssembly< MassAssemblyDim > * | mass_assembly_ |
general assembly objects, hold assembly objects of appropriate dimension More... | |
GenericAssembly< StiffnessAssemblyDim > * | stiffness_assembly_ |
GenericAssembly< SourcesAssemblyDim > * | sources_assembly_ |
GenericAssembly< BdrConditionAssemblyDim > * | bdr_cond_assembly_ |
GenericAssemblyBase * | init_assembly_ |
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 134 of file transport_dg.hh.
using TransportDG< Model >::BdrConditionAssemblyDim = BdrConditionAssemblyDG<dim, Model> |
Definition at line 141 of file transport_dg.hh.
typedef std::vector<std::shared_ptr<FieldFE< 3, FieldValue<3>::Scalar> > > TransportDG< Model >::FieldFEScalarVec |
Definition at line 145 of file transport_dg.hh.
using TransportDG< Model >::InitConditionAssemblyDim = InitConditionAssemblyDG<dim, Model> |
Definition at line 142 of file transport_dg.hh.
using TransportDG< Model >::InitProjectionAssemblyDim = InitProjectionAssemblyDG<dim, Model> |
Definition at line 143 of file transport_dg.hh.
using TransportDG< Model >::MassAssemblyDim = MassAssemblyDG<dim, Model> |
Definition at line 138 of file transport_dg.hh.
using TransportDG< Model >::SourcesAssemblyDim = SourcesAssemblyDG<dim, Model> |
Definition at line 140 of file transport_dg.hh.
using TransportDG< Model >::StiffnessAssemblyDim = StiffnessAssemblyDG<dim, Model> |
Definition at line 139 of file transport_dg.hh.
enum TransportDG::DGVariant |
Enumerator | |
---|---|
non_symmetric | |
incomplete | |
symmetric |
Definition at line 244 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 207 of file transport_dg.cc.
|
override |
Destructor.
Definition at line 357 of file transport_dg.cc.
|
inline |
Access to balance object of Model.
Definition at line 313 of file transport_dg.hh.
void TransportDG< Model >::calculate_cumulative_balance |
Definition at line 650 of file transport_dg.cc.
void TransportDG< Model >::compute_p0_interpolation |
Compute P0 interpolation of the solution (used in reaction term).
Definition at line 598 of file transport_dg.cc.
|
inline |
Definition at line 319 of file transport_dg.hh.
|
inline |
Definition at line 317 of file transport_dg.hh.
|
inline |
Definition at line 277 of file transport_dg.hh.
|
inline |
Return PETSc vector with solution for sbi-th component.
Definition at line 300 of file transport_dg.hh.
|
static |
Input type for the DG variant selection.
Definition at line 53 of file transport_dg.cc.
|
static |
Declare input record type for the equation TransportDG.
Definition at line 79 of file transport_dg.cc.
|
inline |
Getter for P0 interpolation by FieldFE.
Definition at line 304 of file transport_dg.hh.
|
override |
Definition at line 247 of file transport_dg.cc.
void TransportDG< Model >::output_data |
Postprocesses the solution and writes to output file.
Definition at line 625 of file transport_dg.cc.
|
private |
|
private |
Definition at line 432 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 the initial condition.
This method just calls AssemblyDG::prepare_initial_condition() for each elements.
Definition at line 671 of file transport_dg.cc.
void TransportDG< Model >::update_after_reactions | ( | bool | solution_changed | ) |
Definition at line 697 of file transport_dg.cc.
|
override |
Computes the solution in one time instant.
Definition at line 463 of file transport_dg.cc.
|
override |
Initialize solution in the zero time.
Definition at line 399 of file transport_dg.cc.
|
private |
Indicates whether matrices have been preallocated.
Definition at line 415 of file transport_dg.hh.
|
private |
Definition at line 424 of file transport_dg.hh.
|
private |
Data for model parameters.
Definition at line 362 of file transport_dg.hh.
|
private |
Fields for model parameters.
Definition at line 359 of file transport_dg.hh.
|
private |
Definition at line 425 of file transport_dg.hh.
|
private |
Definition at line 417 of file transport_dg.hh.
|
private |
Array for storing the output solution data.
Record with input specification.
Definition at line 392 of file transport_dg.hh.
|
private |
general assembly objects, hold assembly objects of appropriate dimension
Definition at line 421 of file transport_dg.hh.
|
private |
The mass matrix.
Definition at line 378 of file transport_dg.hh.
|
private |
Mass from previous time instant (necessary when coefficients of mass matrix change in time).
Definition at line 381 of file transport_dg.hh.
|
private |
Definition at line 394 of file transport_dg.hh.
|
staticprivate |
Registrar of class to factory.
Definition at line 323 of file transport_dg.hh.
|
private |
Temporary values of increments due to retardation (e.g. sorption)
Definition at line 404 of file transport_dg.hh.
|
private |
Definition at line 404 of file transport_dg.hh.
|
private |
Vector of right hand side.
Definition at line 372 of file transport_dg.hh.
|
private |
Definition at line 423 of file transport_dg.hh.
|
private |
Definition at line 422 of file transport_dg.hh.
|
private |
The stiffness matrix.
Definition at line 375 of file transport_dg.hh.