19 #ifndef TRANSPORT_DG_HH_
20 #define TRANSPORT_DG_HH_
33 template<
unsigned int dim,
unsigned int spacedim>
class FEValuesBase;
34 template<
unsigned int dim,
unsigned int spacedim>
class FiniteElement;
35 template<
unsigned int dim,
unsigned int spacedim>
class Mapping;
51 template<
unsigned int dim>
54 template<
unsigned int dim>
57 template<
unsigned int dim>
60 template<
unsigned int dim>
129 template<
class Model>
134 class EqData :
public Model::ModelEqData {
227 inline typename Model::ModelEqData &
data() {
return data_; }
244 template<
unsigned int dim>
259 template<
unsigned int dim>
272 template<
unsigned int dim>
278 template<
unsigned int dim>
284 template<
unsigned int dim>
290 template<
unsigned int dim>
305 template<
unsigned int dim>
316 template<
unsigned int dim>
331 double Dm,
double alphaL,
double alphaT,
double porosity,
363 double &transport_flux);
396 template<
unsigned int dim>
Input::Record input_rec
Record with input specification.
void assemble_fluxes_element_element()
Assembles the fluxes between elements of the same dimension.
void set_sources()
Assembles the right hand side due to volume sources.
void set_boundary_conditions()
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
FiniteElement< dim, 3 > * fe()
vector< double > mm_coef
Mass matrix coefficients.
void assemble_fluxes_element_side()
Assembles the fluxes between elements of different dimensions.
DOFHandlerMultiDim * dh_
Object for distribution of dofs.
void calculate_concentration_matrix()
Transport with dispersion implemented using discontinuous Galerkin method.
Field< 3, FieldValue< 3 >::Integer > region_id
void calculate_dispersivity_tensor(arma::mat33 &K, const arma::vec3 &velocity, double Dm, double alphaL, double alphaT, double porosity, double cross_cut)
Calculates the dispersivity (diffusivity) tensor from the velocity field.
void assemble_mass_matrix()
Assembles the mass matrix.
FiniteElement< 2, 3 > * fe_rt2_
Wrappers for linear systems based on MPIAIJ and MATIS format.
Class template representing a field with values dependent on: point, element, and region...
void update_solution() override
Computes the solution in one time instant.
FiniteElement< 2, 3 > * fe2_
int dg_variant
DG variant ((non-)symmetric/incomplete.
void prepare_initial_condition()
Assembles the auxiliary linear system to calculate the initial solution as L^2-projection of the pres...
void set_initial_condition()
Sets the initial condition.
FiniteElement< 3, 3 > * fe3_
Mapping< 0, 3 > * map0_
Auxiliary mappings of reference elements.
vector< vector< arma::vec3 > > ad_coef
Advection coefficients.
Mat * mass_matrix
The mass matrix.
void assemble_stiffness_matrix()
Assembles the stiffness matrix.
vector< vector< vector< arma::mat33 > > > dif_coef_edg
Diffusion coefficients on edges.
const Vec & get_solution(unsigned int sbi)
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
static const Input::Type::Selection & get_dg_variant_selection_input_type()
Input type for the DG variant selection.
EquationOutput output_fields
Base class for quadrature rules on simplices in arbitrary dimensions.
void calculate_cumulative_balance()
void assemble_fluxes_boundary()
Assembles the fluxes on the boundary.
TransportDG(Mesh &init_mesh, const Input::Record in_rec)
Constructor.
vector< vector< double > > ret_coef
Retardation coefficient due to sorption.
vector< double > sorption_sources
Temporary values of source increments.
void initialize() override
vector< vector< vector< arma::vec3 > > > ad_coef_edg
Advection coefficients on edges.
void output_data()
Postprocesses the solution and writes to output file.
double ** get_concentration_matrix()
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
vector< vector< arma::mat33 > > dif_coef
Diffusion coefficients.
FEObjects * feo
Finite element objects.
static const int registrar
Registrar of class to factory.
FiniteElement< 1, 3 > * fe_rt1_
Finite elements for the water velocity field.
Abstract class for the mapping between reference and actual cell.
FiniteElement< 1, 3 > * fe1_
Finite elements for the solution of the advection-diffusion equation.
const Vec & get_solution()
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Mat * stiffness_matrix
The stiffness matrix.
bool evaluate_time_constraint(double &time_constraint)
FiniteElement< dim, 3 > * fe_rt()
unsigned int dg_order
Polynomial order of finite elements.
void output_vector_gather()
The class for outputting data during time.
Vec * rhs
Vector of right hand side.
DOFHandlerMultiDim * dh()
double ** solution_elem_
Element averages of solution (the array is passed to reactions in operator splitting).
Model::ModelEqData & data()
FiniteElement< 3, 3 > * fe_rt3_
EqData data_
Field data for model parameters.
bool allocation_done
Indicates whether matrices have been preallocated.
std::vector< std::vector< double > > gamma
Penalty parameters.
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
void update_after_reactions(bool solution_changed)
Vec * mass_vec
Mass from previous time instant (necessary when coefficients of mass matrix change in time)...
Abstract linear system class.
vector< Vec > output_vec
Array for storing the output solution data.
void calculate_instant_balance()
void assemble_volume_integrals()
Assembles the volume integrals into the stiffness matrix.
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...
Base class for FEValues and FESideValues.
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.
Mapping< dim, 3 > * mapping()
~TransportDG()
Destructor.
LinSys ** ls
Linear algebra system for the transport equation.
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)
Sets up some parameters of the DG method for two sides of an edge.
LinSys ** ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
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.
Quadrature< 0 > * q0_
Quadratures used in assembling methods.
void get_par_info(int *&el_4_loc, Distribution *&el_ds)
void zero_time_step() override
Initialize solution in the zero time.
FEObjects(Mesh *mesh_, unsigned int fe_order)