19 #ifndef TRANSPORT_DG_HH_ 20 #define TRANSPORT_DG_HH_ 57 template<
unsigned int dim,
class Model>
class AssemblyDG;
64 template<
unsigned int dim,
unsigned int spacedim>
class FEValuesBase;
66 template<
unsigned int dim,
unsigned int spacedim>
class Mapping;
68 namespace Input {
namespace Type {
class Selection; } }
92 std::shared_ptr<DOFHandlerMultiDim> dh;
131 template<
class Model>
158 class EqData :
public Model::ModelEqData {
176 void set_DG_parameters_boundary(
Side side,
230 std::shared_ptr<DOFHandlerMultiDim>
dh_;
233 std::shared_ptr<DOFHandlerMultiDim>
dh_p0;
269 void zero_time_step()
override;
277 void update_solution()
override;
289 void initialize()
override;
291 void calculate_cumulative_balance();
295 {
return eq_data_->ls[sbi]->get_solution(); }
299 {
return eq_data_->conc_fe;}
302 void compute_p0_interpolation();
304 void update_after_reactions(
bool solution_changed);
311 inline std::shared_ptr<Balance>
balance()
const {
312 return Model::balance_;
315 inline typename Model::ModelEqFields &
eq_fields() {
return *eq_fields_; }
317 inline typename Model::ModelEqData &
eq_data() {
return *eq_data_; }
345 void set_initial_condition();
349 void output_region_statistics();
Input::Record input_rec
Record with input specification.
GenericAssembly< SourcesAssemblyDim > * sources_assembly_
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > FieldFEScalarVec
std::vector< std::vector< double > > gamma
Penalty parameters.
Transport with dispersion implemented using discontinuous Galerkin method.
GenericAssembly< StiffnessAssemblyDim > * stiffness_assembly_
unsigned int dg_order
Polynomial order of finite elements.
GenericAssembly< MassAssemblyDim > * mass_assembly_
general assembly objects, hold assembly objects of appropriate dimension
Wrappers for linear systems based on MPIAIJ and MATIS format.
Class template representing a field with values dependent on: point, element, and region...
vector< VectorMPI > output_vec
Array for storing the output solution data.
Field< 3, FieldValue< 3 >::Scalar > region_id
Vec get_component_vec(unsigned int sbi)
Return PETSc vector with solution for sbi-th component.
FieldFEScalarVec & get_p0_interpolation()
Getter for P0 interpolation by FieldFE.
Directing class of FieldValueCache.
int dg_variant
DG variant ((non-)symmetric/incomplete.
std::shared_ptr< EqFields > eq_fields_
Fields for model parameters.
std::shared_ptr< Balance > balance_
Discontinuous Galerkin method for equation of transport with dispersion.
Basic time management functionality for unsteady (and steady) solvers (class Equation).
vector< double > ret_sources_prev
std::shared_ptr< DOFHandlerMultiDim > dh_p0
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
std::vector< Mat > mass_matrix
The mass matrix.
Base class for quadrature rules on simplices in arbitrary dimensions.
Model::ModelEqFields & eq_fields()
bool evaluate_time_constraint(double &)
unsigned int max_edg_sides
Maximal number of edge sides (evaluate from dim 1,2,3)
vector< vector< arma::mat33 > > dif_coef
Diffusion coefficients.
static const int registrar
Registrar of class to factory.
std::vector< Vec > ret_vec
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
LinSys ** ls
Linear algebra system for the transport equation.
GenericAssembly< BdrConditionAssemblyDim > * bdr_cond_assembly_
Provides the numbering of the finite element degrees of freedom on the computational mesh...
std::vector< Vec > rhs
Vector of right hand side.
Field< 3, FieldValue< 3 >::Scalar > subdomain
The class for outputting data during time.
GenericAssembly< InitConditionAssemblyDim > * init_cond_assembly_
unsigned int IntDim
A dimension index type.
std::vector< Mat > stiffness_matrix
The stiffness matrix.
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
bool allocation_done
Indicates whether matrices have been preallocated.
std::shared_ptr< EqData > eq_data_
Data for model parameters.
void set_time_governor(TimeGovernor *time)
Generic class of assemblation.
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
std::shared_ptr< DOFHandlerMultiDim > dh_
Object for distribution of dofs.
std::shared_ptr< Balance > balance() const
Access to balance object of Model.
Discontinuous Galerkin method for equation of transport with dispersion.
std::vector< Vec > mass_vec
Mass from previous time instant (necessary when coefficients of mass matrix change in time)...
Class for representation of a vector of fields of the same physical quantity.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
#define ASSERT_PTR_DBG(ptr)
Definition of assert macro checking non-null pointer (PTR) only for debug mode.
Model::ModelEqData & eq_data()
LinSys ** ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
EquationOutput output_fields