Flow123d  jenkins-Flow123d-linux-release-multijob-282
Classes | Public Types | Public Member Functions | Static Public Attributes | Private Member Functions | 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 }
 

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...
 
void update_solution () override
 Computes the solution in one time instant. More...
 
virtual void set_velocity_field (const MH_DofHandler &dh)
 Updates the velocity field which determines some coefficients of the transport equation. More...
 
void output_data ()
 Postprocesses the solution and writes to output file. More...
 
virtual EqDataget_data ()
 Getter for field data. More...
 
 ~TransportDG ()
 Destructor. More...
 
template<unsigned int dim>
void calculate_velocity (const typename DOFHandlerBase::CellIterator &cell, vector< arma::vec3 > &velocity, FEValuesBase< dim, 3 > &fv)
 
- Public Member Functions inherited from TransportBase
 TransportBase (Mesh &mesh, const Input::Record in_rec)
 
virtual ~TransportBase ()
 
unsigned int n_substances () override
 Returns number of trnasported substances. More...
 
SubstanceListsubstances () override
 Returns reference to the vector of substnace names. More...
 
virtual void set_concentration_vector (Vec &vec)
 
- Public Member Functions inherited from AdvectionProcessBase
 AdvectionProcessBase (Mesh &mesh, const Input::Record in_rec)
 
- Public Member Functions inherited from EquationBase
 EquationBase ()
 
 EquationBase (Mesh &mesh, const Input::Record in_rec)
 
virtual ~EquationBase ()
 
virtual void initialize ()
 Initialize fields. More...
 
virtual void choose_next_time ()
 
virtual void set_time_upper_constraint (double dt)
 
virtual void set_time_lower_constraint (double dt)
 
TimeGovernor const & time ()
 
virtual void set_time_governor (TimeGovernor &time)
 
double planned_time ()
 
double solved_time ()
 
Meshmesh ()
 
TimeMark::Type mark_type ()
 
FieldSetdata ()
 
virtual void get_solution_vector (double *&vector, unsigned int &size)
 
virtual void get_parallel_solution_vector (Vec &vector)
 

Static Public Attributes

static Input::Type::Record input_type
 Declare input record type for the equation TransportDG. More...
 
static Input::Type::Selection dg_variant_selection_input_type
 Input type for the DG variant selection. More...
 
- Static Public Attributes inherited from TransportBase
static Input::Type::Record input_type_output_record
 
- Static Public Attributes inherited from AdvectionProcessBase
static Input::Type::AbstractRecord input_type
 Common specification of the input record for secondary equations. More...
 

Private Member Functions

Model::ModelEqData & data ()
 
void output_vector_gather ()
 
void preallocate ()
 
void assemble_mass_matrix ()
 Assembles the mass matrix. More...
 
template<unsigned int dim>
void assemble_mass_matrix ()
 Assembles the mass matrix for the given dimension. More...
 
void assemble_stiffness_matrix ()
 Assembles the stiffness matrix. More...
 
template<unsigned int dim>
void assemble_volume_integrals ()
 Assembles the volume integrals into the stiffness matrix. More...
 
void set_sources ()
 Assembles the right hand side due to volume sources. More...
 
template<unsigned int dim>
void set_sources ()
 Assembles the right hand side vector due to volume sources. More...
 
template<unsigned int dim>
void assemble_fluxes_boundary ()
 Assembles the fluxes on the boundary. More...
 
template<unsigned int dim>
void assemble_fluxes_element_element ()
 Assembles the fluxes between elements of the same dimension. More...
 
template<unsigned int dim>
void assemble_fluxes_element_side ()
 Assembles the fluxes between elements of different dimensions. More...
 
void set_boundary_conditions ()
 Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions. More...
 
template<unsigned int dim>
void set_boundary_conditions ()
 Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions for a given space dimension. More...
 
template<unsigned int dim>
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. More...
 
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. More...
 
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. More...
 
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. More...
 
void set_initial_condition ()
 Sets the initial condition. More...
 
template<unsigned int dim>
void prepare_initial_condition ()
 Assembles the auxiliary linear system to calculate the initial solution as L^2-projection of the prescribed initial condition. More...
 
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. More...
 
template<unsigned int dim>
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 of specific dimension. More...
 
void calc_elem_sources (vector< vector< double > > &mass, vector< vector< double > > &src_balance)
 Calculates volume sources for each region. More...
 
template<unsigned int dim>
void calc_elem_sources (vector< vector< double > > &mass, vector< vector< double > > &src_balance)
 Calculates volume sources for each region of specific dimension. More...
 

Private Attributes

Physical parameters
EqData data_
 Field data for model parameters. More...
 
Parameters of the numerical method
FEObjectsfeo
 Finite element objects. More...
 
std::vector< std::vector
< double > > 
gamma
 Penalty parameters. More...
 
int dg_variant
 DG variant ((non-)symmetric/incomplete. More...
 
unsigned int dg_order
 Polynomial order of finite elements. More...
 
Solution of algebraic system
Vec * rhs
 Vector of right hand side. More...
 
Mat * stiffness_matrix
 The stiffness matrix. More...
 
Mat * mass_matrix
 The mass matrix. More...
 
Vec * mass_vec
 Mass from previous time instant (necessary when coefficients of mass matrix change in time). More...
 
LinSys ** ls
 Linear algebra system for the transport equation. More...
 
LinSys ** ls_dt
 Linear algebra system for the time derivative (actually it is used only for handling the matrix structures). More...
 
Output to file
vector< double * > output_solution
 Array for storing the output solution data. More...
 
vector< Vec > output_vec
 Vector of solution data. More...
 
Input::Record output_rec
 Record with output specification. More...
 
OutputTimeoutput_stream
 
Auxiliary fields used during assembly
vector< double > mm_coef
 Mass matrix coefficients. More...
 
vector< vector< double > > ret_coef
 Retardation coefficient due to sorption. More...
 
vector< double > sorption_sources
 Temporary values of source increments. More...
 
vector< vector< arma::vec3 > > ad_coef
 Advection coefficients. More...
 
vector< vector< arma::mat33 > > dif_coef
 Diffusion coefficients. More...
 
vector< vector< vector
< arma::vec3 > > > 
ad_coef_edg
 Advection coefficients on edges. More...
 
vector< vector< vector
< arma::mat33 > > > 
dif_coef_edg
 Diffusion coefficients on edges. More...
 
vector< unsigned int > subst_idx
 List of indices used to call balance methods for a set of quantities. More...
 
Other
bool allocation_done
 Indicates whether matrices have been preallocated. More...
 

Additional Inherited Members

- Protected Member Functions inherited from TransportBase
const RegionDBregion_db ()
 Returns the region database. More...
 
- Protected Attributes inherited from TransportBase
unsigned int n_subst_
 Number of transported substances. More...
 
SubstanceList substances_
 Transported substances. More...
 
const MH_DofHandlermh_dh
 
boost::shared_ptr< Balancebalance_
 (new) object for calculation and writing the mass balance to file. More...
 
- Protected Attributes inherited from EquationBase
bool equation_empty_
 flag is true if only default constructor was called More...
 
Meshmesh_
 
TimeGovernortime_
 
Input::Record input_record_
 
FieldSeteq_data_
 

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 137 of file transport_dg.hh.

Member Enumeration Documentation

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

Definition at line 170 of file transport_dg.hh.

Constructor & Destructor Documentation

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

Constructor.

Parameters
init_meshcomputational mesh
in_recinput record

Definition at line 249 of file transport_dg.cc.

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

Destructor.

Definition at line 387 of file transport_dg.cc.

Member Function Documentation

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::assemble_fluxes_boundary ( )
private

Assembles the fluxes on the boundary.

Definition at line 1068 of file transport_dg.cc.

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::assemble_fluxes_element_element ( )
private

Assembles the fluxes between elements of the same dimension.

Definition at line 937 of file transport_dg.cc.

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::assemble_fluxes_element_side ( )
private

Assembles the fluxes between elements of different dimensions.

Definition at line 1161 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::assemble_mass_matrix ( )
private

Assembles the mass matrix.

The routine just calls templated method assemble_mass_matrix() for each space dimension.

Definition at line 705 of file transport_dg.cc.

template<class Model>
template<unsigned int dim>
void TransportDG< Model >::assemble_mass_matrix ( )
private

Assembles the mass matrix for the given dimension.

template<class Model >
void TransportDG< Model >::assemble_stiffness_matrix ( )
private

Assembles the stiffness matrix.

This routine just calls assemble_volume_integrals(), assemble_fluxes_boundary(), assemble_fluxes_element_element() and assemble_fluxes_element_side() for each space dimension.

Definition at line 773 of file transport_dg.cc.

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::assemble_volume_integrals ( )
private

Assembles the volume integrals into the stiffness matrix.

Definition at line 806 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::calc_elem_sources ( vector< vector< double > > &  mass,
vector< vector< double > > &  src_balance 
)
private

Calculates volume sources for each region.

This method actually calls calc_elem_sources<dim>() for each space dimension.

Parameters
massVector of substance mass per region.
src_balanceVector of sources per region.

Definition at line 1773 of file transport_dg.cc.

template<class Model>
template<unsigned int dim>
void TransportDG< Model >::calc_elem_sources ( vector< vector< double > > &  mass,
vector< vector< double > > &  src_balance 
)
private

Calculates volume sources for each region of specific dimension.

Parameters
massVector of substance mass per region.
src_balanceVector of sources per region.
template<class Model >
void TransportDG< Model >::calc_fluxes ( vector< vector< double > > &  bcd_balance,
vector< vector< double > > &  bcd_plus_balance,
vector< vector< double > > &  bcd_minus_balance 
)
private

Calculates flux through boundary of each region.

This actually calls calc_fluxes<dim>() for each space dimension.

Parameters
bcd_balanceTotal fluxes.
bcd_plus_balanceIncoming fluxes.
bcd_minus_balanceOutgoing fluxes.

Definition at line 1674 of file transport_dg.cc.

template<class Model>
template<unsigned int dim>
void TransportDG< Model >::calc_fluxes ( vector< vector< double > > &  bcd_balance,
vector< vector< double > > &  bcd_plus_balance,
vector< vector< double > > &  bcd_minus_balance 
)
private

Calculates flux through boundary of each region of specific dimension.

Parameters
bcd_balanceTotal fluxes.
bcd_plus_balanceIncoming fluxes.
bcd_minus_balanceOutgoing fluxes.
template<class Model>
void TransportDG< Model >::calculate_dispersivity_tensor ( arma::mat33 K,
const arma::vec3 velocity,
double  Dm,
double  alphaL,
double  alphaT,
double  porosity,
double  cross_cut 
)
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.
template<class Model>
template<unsigned int dim>
void TransportDG< Model >::calculate_velocity ( const ElementFullIter cell,
std::vector< arma::vec3 > &  velocity,
FEValuesBase< dim, 3 > &  fv 
)
private

Calculates the velocity field on a given dim dimensional cell.

Parameters
cellThe cell.
velocityThe computed velocity field (at quadrature points).
fvThe FEValues class providing the quadrature points and the shape functions for velocity.
template<class Model>
template<unsigned int dim>
void TransportDG< Model >::calculate_velocity ( const typename DOFHandlerBase::CellIterator cell,
vector< arma::vec3 > &  velocity,
FEValuesBase< dim, 3 > &  fv 
)

Definition at line 1434 of file transport_dg.cc.

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

Definition at line 235 of file transport_dg.hh.

template<class Model>
virtual EqData* TransportDG< Model >::get_data ( )
inlinevirtual

Getter for field data.

Definition at line 226 of file transport_dg.hh.

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

Postprocesses the solution and writes to output file.

Reimplemented from EquationBase.

Definition at line 685 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::output_vector_gather ( void  )
private

Definition at line 423 of file transport_dg.cc.

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

Definition at line 484 of file transport_dg.cc.

template<class Model >
template<unsigned int dim>
void TransportDG< Model >::prepare_initial_condition ( )
private

Assembles the auxiliary linear system to calculate the initial solution as L^2-projection of the prescribed initial condition.

Definition at line 1623 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::set_boundary_conditions ( )
private

Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.

The routine just calls templated method set_boundary_condition() for each space dimension.

Definition at line 1275 of file transport_dg.cc.

template<class Model>
template<unsigned int dim>
void TransportDG< Model >::set_boundary_conditions ( )
private

Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions for a given space dimension.

template<class Model >
void TransportDG< Model >::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 
)
private

Sets up parameters of the DG method on a given boundary edge.

Assumption is that the edge consists of only 1 side.

Parameters
sideThe boundary side.
K_sizeSize of vector of tensors K.
KDispersivity tensor.
ad_vectorAdvection vector.
normal_vectorNormal vector (assumed constant along the edge).
alphaPenalty parameter that influences the continuity of the solution (large value=more continuity).
gammaComputed penalty parameters.

Definition at line 1563 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::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 
)
private

Sets up some parameters of the DG method for two sides of an edge.

Parameters
edgThe edge.
s1Side 1.
s2Side 2.
K_sizeSize of vector of tensors K.
K1Dispersivity tensors on side s1 (in quadrature points).
K2Dispersivity tensors on side s2 (in quadrature points).
normal_vectorNormal vector to side 0 of the neighbour (assumed constant along the side).
alpha1,alpha2Penalty parameter that influences the continuity of the solution (large value=more continuity).
gammaComputed penalty parameters.
omegaComputed weights.
transport_fluxComputed flux from side s1 to side s2.

Definition at line 1461 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::set_initial_condition ( )
private

Sets the initial condition.

Definition at line 1598 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::set_sources ( )
private

Assembles the right hand side due to volume sources.

This method just calls set_sources() for each space dimension.

Definition at line 860 of file transport_dg.cc.

template<class Model>
template<unsigned int dim>
void TransportDG< Model >::set_sources ( )
private

Assembles the right hand side vector due to volume sources.

template<class Model >
void TransportDG< Model >::set_velocity_field ( const MH_DofHandler dh)
virtual

Updates the velocity field which determines some coefficients of the transport equation.

@param dh mixed hybrid dof handler

(So far it does not work since the flow module returns a vector of zeros.)

Parameters
velocity_vectorInput array of velocity values.

Reimplemented from TransportBase.

Definition at line 673 of file transport_dg.cc.

template<class Model >
void TransportDG< Model >::update_solution ( void  )
overridevirtual

Computes the solution in one time instant.

Reimplemented from EquationBase.

Definition at line 508 of file transport_dg.cc.

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

Initialize solution in the zero time.

Reimplemented from EquationBase.

Definition at line 443 of file transport_dg.cc.

Member Data Documentation

template<class Model>
vector<vector<arma::vec3> > TransportDG< Model >::ad_coef
private

Advection coefficients.

Definition at line 525 of file transport_dg.hh.

template<class Model>
vector<vector<vector<arma::vec3> > > TransportDG< Model >::ad_coef_edg
private

Advection coefficients on edges.

Definition at line 529 of file transport_dg.hh.

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

Indicates whether matrices have been preallocated.

Definition at line 544 of file transport_dg.hh.

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

Field data for model parameters.

Definition at line 449 of file transport_dg.hh.

template<class Model>
unsigned int TransportDG< Model >::dg_order
private

Polynomial order of finite elements.

Definition at line 467 of file transport_dg.hh.

template<class Model>
int TransportDG< Model >::dg_variant
private

DG variant ((non-)symmetric/incomplete.

Definition at line 464 of file transport_dg.hh.

template<class Model>
Selection TransportDG< Model >::dg_variant_selection_input_type
static
Initial value:
= Selection("DG_variant", "Type of penalty term.")
.add_value(non_symmetric, "non-symmetric", "non-symmetric weighted interior penalty DG method")
.add_value(incomplete, "incomplete", "incomplete weighted interior penalty DG method")
.add_value(symmetric, "symmetric", "symmetric weighted interior penalty DG method")

Input type for the DG variant selection.

Definition at line 196 of file transport_dg.hh.

template<class Model>
vector<vector<arma::mat33> > TransportDG< Model >::dif_coef
private

Diffusion coefficients.

Definition at line 527 of file transport_dg.hh.

template<class Model>
vector<vector<vector<arma::mat33> > > TransportDG< Model >::dif_coef_edg
private

Diffusion coefficients on edges.

Definition at line 531 of file transport_dg.hh.

template<class Model>
FEObjects* TransportDG< Model >::feo
private

Finite element objects.

Definition at line 458 of file transport_dg.hh.

template<class Model>
std::vector<std::vector<double> > TransportDG< Model >::gamma
private

Penalty parameters.

Definition at line 461 of file transport_dg.hh.

template<class Model>
Record TransportDG< Model >::input_type
static
Initial value:
= Model::get_input_type("DG", "DG solver")
.declare_key("solver", LinSys_PETSC::input_type, Default::obligatory(),
"Linear solver for MH problem.")
.declare_key("input_fields", Array(TransportDG<Model>::EqData().make_field_descriptor_type(std::string(Model::ModelEqData::name()) + "_DG")), IT::Default::obligatory(), "")
.declare_key("dg_variant", TransportDG<Model>::dg_variant_selection_input_type, Default("non-symmetric"),
"Variant of interior penalty discontinuous Galerkin method.")
.declare_key("dg_order", Integer(0,3), Default("1"),
"Polynomial order for finite element in DG method (order 0 is suitable if there is no diffusion/dispersion).")
.declare_key("output_fields", Array(EqData::output_selection),
Default(Model::ModelEqData::default_output_field()),
"List of fields to write to output file.")

Declare input record type for the equation TransportDG.

Definition at line 191 of file transport_dg.hh.

template<class Model>
LinSys** TransportDG< Model >::ls
private

Linear algebra system for the transport equation.

Definition at line 489 of file transport_dg.hh.

template<class Model>
LinSys** TransportDG< Model >::ls_dt
private

Linear algebra system for the time derivative (actually it is used only for handling the matrix structures).

Definition at line 492 of file transport_dg.hh.

template<class Model>
Mat* TransportDG< Model >::mass_matrix
private

The mass matrix.

Definition at line 483 of file transport_dg.hh.

template<class Model>
Vec* TransportDG< Model >::mass_vec
private

Mass from previous time instant (necessary when coefficients of mass matrix change in time).

Definition at line 486 of file transport_dg.hh.

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

Mass matrix coefficients.

Definition at line 519 of file transport_dg.hh.

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

Record with output specification.

Definition at line 507 of file transport_dg.hh.

template<class Model>
vector<double*> TransportDG< Model >::output_solution
private

Array for storing the output solution data.

Definition at line 501 of file transport_dg.hh.

template<class Model>
OutputTime* TransportDG< Model >::output_stream
private

Definition at line 509 of file transport_dg.hh.

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

Vector of solution data.

Definition at line 504 of file transport_dg.hh.

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

Retardation coefficient due to sorption.

Definition at line 521 of file transport_dg.hh.

template<class Model>
Vec* TransportDG< Model >::rhs
private

Vector of right hand side.

Definition at line 477 of file transport_dg.hh.

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

Temporary values of source increments.

Definition at line 523 of file transport_dg.hh.

template<class Model>
Mat* TransportDG< Model >::stiffness_matrix
private

The stiffness matrix.

Definition at line 480 of file transport_dg.hh.

template<class Model>
vector<unsigned int> TransportDG< Model >::subst_idx
private

List of indices used to call balance methods for a set of quantities.

Definition at line 533 of file transport_dg.hh.


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