30 #ifndef TRANSPORT_DG_HH_
31 #define TRANSPORT_DG_HH_
40 template<
unsigned int dim,
unsigned int spacedim>
class FEValuesBase;
41 template<
unsigned int dim,
unsigned int spacedim>
class FiniteElement;
42 template<
unsigned int dim,
unsigned int spacedim>
class Mapping;
58 template<
unsigned int dim>
61 template<
unsigned int dim>
64 template<
unsigned int dim>
67 template<
unsigned int dim>
136 template<
class Model>
141 class EqData :
public Model::ModelEqData {
234 inline typename Model::ModelEqData &
data() {
return data_; }
249 template<
unsigned int dim>
264 template<
unsigned int dim>
277 template<
unsigned int dim>
283 template<
unsigned int dim>
289 template<
unsigned int dim>
295 template<
unsigned int dim>
310 template<
unsigned int dim>
321 template<
unsigned int dim>
336 double Dm,
double alphaL,
double alphaT,
double porosity,
368 double &transport_flux);
401 template<
unsigned int dim>
420 template<
unsigned int dim>
437 template<
unsigned int dim>
static Input::Type::Selection output_selection
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.
Input::Record output_rec
Record with output specification.
Transport with dispersion implemented using discontinuous Galerkin method.
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_
OutputTime * output_stream
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.
void assemble_stiffness_matrix()
Assembles the stiffness matrix.
vector< vector< vector< arma::mat33 > > > dif_coef_edg
Diffusion coefficients on edges.
Specification of transport model interface.
Base class for quadrature rules on simplices in arbitrary dimensions.
void assemble_fluxes_boundary()
Assembles the fluxes on the boundary.
TransportDG(Mesh &init_mesh, const Input::Record &in_rec)
Constructor.
vector< double * > output_solution
Array for storing the output solution data.
vector< vector< vector< arma::vec3 > > > ad_coef_edg
Advection coefficients on edges.
void output_data()
Postprocesses the solution and writes to output file.
vector< vector< arma::mat33 > > dif_coef
Diffusion coefficients.
FEObjects * feo
Finite element objects.
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.
Mat mass_matrix
The mass matrix.
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Mat * stiffness_matrix
The stiffness matrix.
Field< 3, FieldValue< 3 >::Vector > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
FiniteElement< dim, 3 > * fe_rt()
unsigned int dg_order
Polynomial order of finite elements.
void output_vector_gather()
The class for outputing data during time.
Vec * rhs
Vector of right hand side.
static Input::Type::Record input_type
Declare input record type for the equation TransportDG.
BCField< 3, FieldValue< 3 >::Vector > bc_flux
Flux in Neumann or Robin b.c.
static Input::Type::Selection bc_type_selection
DOFHandlerMultiDim * dh()
LinSys * ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
TimeIntegrationScheme time_scheme()
Returns the time integration scheme of the equation.
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.
void calc_fluxes(vector< vector< double > > &bcd_balance, vector< vector< double > > &bcd_plus_balance, vector< vector< double > > &bcd_minus_balance)
Calculates flux through boundary of each region.
BCField< 3, FieldValue< 3 >::EnumVector > bc_type
Type of boundary condition (see also BC_Type)
Abstract linear system class.
virtual EqData * get_data()
Getter for field data.
void calc_elem_sources(vector< vector< double > > &mass, vector< vector< double > > &src_balance)
Calculates volume sources for each region.
vector< Vec > output_vec
Vector of solution data.
virtual void set_velocity_field(const MH_DofHandler &dh)
Updates the velocity field which determines some coefficients of the transport equation.
void assemble_volume_integrals()
Assembles the volume integrals into the stiffness matrix.
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()
static Input::Type::Selection dg_variant_selection_input_type
Input type for the DG variant selection.
~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.
Field< 3, FieldValue< 3 >::Vector > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
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.
BCField< 3, FieldValue< 3 >::Vector > bc_robin_sigma
Transition coefficient in Robin b.c.
void zero_time_step() override
Initialize solution in the zero time.
FEObjects(Mesh *mesh_, unsigned int fe_order)