Flow123d  JS_before_hm-2158-g8a5161593
transport_dg.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file transport_dg.hh
15  * @brief Discontinuous Galerkin method for equation of transport with dispersion.
16  * @author Jan Stebel
17  */
18 
19 #ifndef TRANSPORT_DG_HH_
20 #define TRANSPORT_DG_HH_
21 
22 #include <math.h> // for fabs
23 #include <string.h> // for memcpy
24 #include <algorithm> // for max
25 
26 #include <string> // for operator<<
27 #include <vector> // for vector
28 #include <armadillo>
29 #include "fem/update_flags.hh" // for operator|
30 #include "fields/field_values.hh" // for FieldValue<>::Scalar
31 #include "fields/field.hh"
32 #include "fields/field_fe.hh"
33 #include "fields/multi_field.hh"
34 #include "la/vector_mpi.hh"
36 #include "la/linsys.hh"
37 #include "input/accessors.hh" // for ExcAccessorForNullStorage
38 #include "input/accessors_impl.hh" // for Record::val
39 #include "input/storage.hh" // for ExcStorageTypeMismatch
40 #include "input/type_base.hh" // for Array
41 #include "input/type_generic.hh" // for Instance
42 #include "input/type_record.hh" // for Record::ExcRecordKeyNo...
43 #include "system/index_types.hh" // for LongIdx
44 #include "mesh/accessors.hh" // for ElementAccessor, SIdeIter
45 #include "mesh/elements.h" // for Element::dim, Element:...
46 #include "mesh/neighbours.h" // for Neighbour::element
47 #include "mpi.h" // for MPI_Comm_rank
48 #include "petscmat.h" // for Mat, MatDestroy
49 #include "petscvec.h" // for Vec, VecDestroy, VecSc...
50 #include "transport/concentration_model.hh" // for ConcentrationTransport...
51 #include "transport/heat_model.hh" // for HeatTransferModel, Hea...
52 #include "tools/mixed.hh"
53 class DiscreteSpace;
54 class Distribution;
55 class OutputTime;
56 class DOFHandlerMultiDim;
58 template<unsigned int dim, class Model> class AssemblyDG;
59 template<unsigned int dim, class Model> class MassAssemblyDG;
60 template<unsigned int dim, class Model> class StiffnessAssemblyDG;
61 template<unsigned int dim, class Model> class SourcesAssemblyDG;
62 template<unsigned int dim, class Model> class BdrConditionAssemblyDG;
63 template<unsigned int dim, class Model> class InitConditionAssemblyDG;
64 template<unsigned int dim, class Model> class InitProjectionAssemblyDG;
65 template< template<IntDim...> class DimAssembly> class GenericAssembly;
66 template<unsigned int dim, unsigned int spacedim> class FEValuesBase;
67 template<unsigned int dim> class FiniteElement;
68 template<unsigned int dim, unsigned int spacedim> class Mapping;
69 class Quadrature;
70 namespace Input { namespace Type { class Selection; } }
71 class ElementCacheMap;
72 
73 /*class FEObjects {
74 public:
75 
76  inline FEObjects(Mesh *mesh_, unsigned int fe_order)
77  {
78  fe = MixedPtr<FE_P_disc>(fe_order);
79  fe_rt = MixedPtr<FE_RT0>();
80  q = MixedPtr<QGauss>(2*fe_order);
81 
82  auto ds = std::make_shared<EqualOrderDiscreteSpace>(mesh_, fe);
83  dh = std::make_shared<DOFHandlerMultiDim>(*mesh_);
84  dh->distribute_dofs(ds_);
85  }
86 
87  ~FEObjects();
88 
89  MixedPtr<FiniteElement> fe;
90  MixedPtr<FiniteElement> fe_rt;
91  MixedPtr<Quadrature> q;
92 
93  /// Object for distribution of dofs.
94  std::shared_ptr<DOFHandlerMultiDim> dh;
95 };*/
96 
97 
98 
99 /**
100  * @brief Transport with dispersion implemented using discontinuous Galerkin method.
101  *
102  * TransportDG implements the discontinuous Galerkin method for the transport and diffusion of substances.
103  * The concentration @f$ c_i ~[kg/m^3]@f$ of the i-th substance is governed by the advection-diffusion equation
104  * @f[
105  * \partial_t c_i + \mathbf v\cdot\nabla c_i - \mathrm{div}(D\nabla c_i) = F \mbox{ in }\Omega^d,
106  * @f]
107  * where @f$\mathbf v@f$ is the fluid velocity and @f$\Omega^d@f$ the @f$d@f$-dimensional domain, respectively.
108  * The hydrodynamic dispersivity tensor @f$\mathbf D ~[m^2/s]@f$ is given by:
109  * @f[
110  * \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).
111  * @f]
112  * The molecular dispersivity @f$D_m~[m^2/s]@f$, as well as the longitudal and transversal dispersivity @f$\alpha_L,~\alpha_T~[m]@f$ are input parameters of the model.
113  *
114  * For lower dimensions @f$d=1,2@f$ the advection-diffusion equation is multiplied by the fracture cross-cut @f$\delta^d~[m^{3-d}]@f$.
115  *
116  * The boundary @f$\partial\Omega^d@f$ is divided into three disjoint parts @f$\Gamma^d_D\cup\Gamma^d_N\cup\Gamma^d_F@f$.
117  * We prescribe the following boundary conditions:
118  * @f{eqnarray*}{
119  * c_i^d &= c_{iD}^d &\mbox{ on }\Gamma^d_D \mbox{ (Dirichlet)},\\
120  * \mathbf D^d\nabla c_i^d\cdot\mathbf n &= 0 &\mbox{ on }\Gamma^d_N \mbox{ (Neumann)},
121  * @f}
122  * The transfer of mass through fractures is described by the transmission conditions on @f$\Gamma^d_F@f$:
123  * @f[
124  * -\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
125  * F^{d-1} = (\sigma + |\mathbf v^d\cdot\mathbf n|)(c_i^d-c_i^{d-1}).
126  * @f]
127  * Here @f$\mathbf n@f$ stands for the unit outward normal vector to @f$\partial\Omega^d@f$.
128  * The coefficient @f$\sigma@f$ determines the transfer of mass through fractures due to diffusion.
129  *
130  * @ingroup transport_mod
131  *
132  */
133 template<class Model>
134 class TransportDG : public Model
135 {
136 public:
137 
138  template<unsigned int dim> using MassAssemblyDim = MassAssemblyDG<dim, Model>;
139  template<unsigned int dim> using StiffnessAssemblyDim = StiffnessAssemblyDG<dim, Model>;
140  template<unsigned int dim> using SourcesAssemblyDim = SourcesAssemblyDG<dim, Model>;
144 
146 
147  class EqFields : public Model::ModelEqFields {
148  public:
149 
150  EqFields();
151 
152  MultiField<3, FieldValue<3>::Scalar> fracture_sigma; ///< Transition parameter for diffusive transfer on fractures (for each substance).
153  MultiField<3, FieldValue<3>::Scalar> dg_penalty; ///< Penalty enforcing inter-element continuity of solution (for each substance).
156 
158  };
159 
160 
161  class EqData : public Model::ModelEqData {
162  public:
163 
164  EqData() {}
165 
166  /**
167  * @brief Sets up parameters of the DG method on a given boundary edge.
168  *
169  * Assumption is that the edge consists of only 1 side.
170  * @param side The boundary side.
171  * @param K_size Size of vector of tensors K.
172  * @param K Dispersivity tensor.
173  * @param ad_vector Advection vector.
174  * @param normal_vector Normal vector (assumed constant along the edge).
175  * @param alpha Penalty parameter that influences the continuity
176  * of the solution (large value=more continuity).
177  * @param gamma Computed penalty parameters.
178  */
180  const int K_size,
181  const std::vector<arma::mat33> &K,
182  const double flux,
183  const arma::vec3 &normal_vector,
184  const double alpha,
185  double &gamma);
186 
187 
188  inline void set_time_governor(TimeGovernor *time) {
189  ASSERT_PTR_DBG(time);
190  this->time_ = time;
191  }
192 
193 
194  /// Compute and return anisotropy of given element
195  double elem_anisotropy(ElementAccessor<3> e) const;
196 
197 
198  /// @name Parameters of the numerical method
199  // @{
200 
201  /// Penalty parameters.
203 
204  /// DG variant ((non-)symmetric/incomplete
206 
207  /// Polynomial order of finite elements.
208  unsigned int dg_order;
209 
210  // @}
211 
212 
213  /// Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
215 
216  /// Linear algebra system for the transport equation.
218 
219  /// Linear algebra system for the time derivative (actually it is used only for handling the matrix structures).
221 
222  /// @name Auxiliary fields used during assembly
223  // @{
224 
225  /// Diffusion coefficients.
227  /// Maximal number of edge sides (evaluate from dim 1,2,3)
228  unsigned int max_edg_sides;
229 
230  // @}
231 
232  /// Object for distribution of dofs.
233  std::shared_ptr<DOFHandlerMultiDim> dh_;
234 
235  /// Vector of solution data.
237 
239  std::shared_ptr<DOFHandlerMultiDim> dh_p0;
241  std::shared_ptr<Balance> balance_;
242  };
243 
244  enum DGVariant {
245  // Non-symmetric weighted interior penalty DG
247 
248  // Incomplete weighted interior penalty DG
250 
251  // Symmetric weighted interior penalty DG
253  };
254 
255  /**
256  * @brief Constructor.
257  * @param init_mesh computational mesh
258  * @param in_rec input record
259  */
260  TransportDG(Mesh &init_mesh, const Input::Record in_rec);
261  /**
262 
263  * @brief Declare input record type for the equation TransportDG.
264  */
265  static const Input::Type::Record & get_input_type();
266 
267  /**
268  * @brief Input type for the DG variant selection.
269  */
271 
272  /**
273  * @brief Initialize solution in the zero time.
274  */
275  void zero_time_step() override;
276 
278  { return false; }
279 
280  /**
281  * @brief Computes the solution in one time instant.
282  */
283  void update_solution() override;
284 
285  /**
286  * @brief Postprocesses the solution and writes to output file.
287  */
288  void output_data();
289 
290  /**
291  * @brief Destructor.
292  */
293  ~TransportDG() override;
294 
295  void initialize() override;
296 
298 
299  /// Return PETSc vector with solution for sbi-th component.
300  Vec get_component_vec(unsigned int sbi)
301  { return eq_data_->ls[sbi]->get_solution(); }
302 
303  /// Getter for P0 interpolation by FieldFE.
305  { return eq_data_->conc_fe;}
306 
307  /// Compute P0 interpolation of the solution (used in reaction term).
309 
310  void update_after_reactions(bool solution_changed);
311 
312  /// Access to balance object of Model
313  inline std::shared_ptr<Balance> balance() const {
314  return Model::balance_;
315  }
316 
317  inline typename Model::ModelEqFields &eq_fields() { return *eq_fields_; }
318 
319  inline typename Model::ModelEqData &eq_data() { return *eq_data_; }
320 
321 private:
322  /// Registrar of class to factory
323  static const int registrar;
324 
325  void preallocate();
326 
327  /**
328  * @brief Calculates the dispersivity (diffusivity) tensor from the velocity field.
329  *
330  * @param K The computed dispersivity tensor.
331  * @param velocity The velocity field (at quadrature points).
332  * @param Dm Molecular diffusivities.
333  * @param alphaL Longitudal dispersivities.
334  * @param alphaT Transversal dispersivities.
335  * @param porosity Porosities.
336  * @param cross_cut Cross-cuts of higher dimension.
337  */
338 // void calculate_dispersivity_tensor(arma::mat33 &K, const arma::vec3 &velocity,
339 // double Dm, double alphaL, double alphaT, double porosity,
340 // double cross_cut);
341 
342  /**
343  * @brief Sets the initial condition.
344  *
345  * This method just calls AssemblyDG::prepare_initial_condition() for each elements.
346  */
347  void set_initial_condition();
348 
349 
350 
352 
353 
354 
355  /// @name Physical parameters
356  // @{
357 
358  /// Fields for model parameters.
359  std::shared_ptr<EqFields> eq_fields_;
360 
361  /// Data for model parameters.
362  std::shared_ptr<EqData> eq_data_;
363 
364  // @}
365 
366 
367 
368  /// @name Solution of algebraic system
369  // @{
370 
371  /// Vector of right hand side.
373 
374  /// The stiffness matrix.
376 
377  /// The mass matrix.
379 
380  /// Mass from previous time instant (necessary when coefficients of mass matrix change in time).
382  // @}
383 
384 
385  /// @name Output to file
386  // @{
387 
388  /// Array for storing the output solution data.
389  //vector<double*> output_solution;
390 
391  /// Record with input specification.
393 
394  ofstream reg_stat_stream;
395 
396 
397  // @}
398 
399 
400  /// @name Auxiliary fields used during assembly
401  // @{
402 
403  /// Temporary values of increments due to retardation (e.g. sorption)
405 
406  // @}
407 
408 
409 
410 
411  /// @name Other
412  // @{
413 
414  /// Indicates whether matrices have been preallocated.
416 
418  // @}
419 
420  /// general assembly objects, hold assembly objects of appropriate dimension
426 };
427 
428 
429 
430 
431 
432 
433 #endif /* TRANSPORT_DG_HH_ */
TransportDG::EqData::dg_order
unsigned int dg_order
Polynomial order of finite elements.
Definition: transport_dg.hh:208
TransportDG::init_projection
bool init_projection
Definition: transport_dg.hh:417
TransportDG::EqData::set_time_governor
void set_time_governor(TimeGovernor *time)
Definition: transport_dg.hh:188
TransportDG::stiffness_assembly_
GenericAssembly< StiffnessAssemblyDim > * stiffness_assembly_
Definition: transport_dg.hh:422
TransportDG::EqFields::dg_penalty
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
Definition: transport_dg.hh:153
TransportDG::get_p0_interpolation
FieldFEScalarVec & get_p0_interpolation()
Getter for P0 interpolation by FieldFE.
Definition: transport_dg.hh:304
TransportDG::FieldFEScalarVec
std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > FieldFEScalarVec
Definition: transport_dg.hh:145
LinSys
Definition: la_linsys_new.hh:169
MassAssemblyDG
Definition: assembly_dg.hh:36
vector_mpi.hh
GenericAssemblyBase
Definition: generic_assembly.hh:52
TransportDG::incomplete
@ incomplete
Definition: transport_dg.hh:249
TransportDG::EqData::conc_fe
FieldFEScalarVec conc_fe
Definition: transport_dg.hh:238
neighbours.h
Input
Abstract linear system class.
Definition: balance.hh:40
TransportDG::EqData::max_edg_sides
unsigned int max_edg_sides
Maximal number of edge sides (evaluate from dim 1,2,3)
Definition: transport_dg.hh:228
string.h
TransportDG::symmetric
@ symmetric
Definition: transport_dg.hh:252
TransportDG::EqData::time_
TimeGovernor * time_
Definition: transport_dg.hh:240
TransportDG::input_rec
Input::Record input_rec
Array for storing the output solution data.
Definition: transport_dg.hh:392
ElementCacheMap
Directing class of FieldValueCache.
Definition: field_value_cache.hh:151
linsys.hh
Wrappers for linear systems based on MPIAIJ and MATIS format.
TransportDG::EqFields::EqFields
EqFields()
Definition: transport_dg.cc:116
std::vector
Definition: doxy_dummy_defs.hh:7
TransportDG::update_solution
void update_solution() override
Computes the solution in one time instant.
Definition: transport_dg.cc:463
ElementAccessor< 3 >
TransportDG::EqData::ls_dt
LinSys ** ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
Definition: transport_dg.hh:220
arma::vec3
Definition: doxy_dummy_defs.hh:17
TransportDG::ret_sources_prev
vector< double > ret_sources_prev
Definition: transport_dg.hh:404
field_fe.hh
type_base.hh
TransportDG::set_initial_condition
void set_initial_condition()
Calculates the dispersivity (diffusivity) tensor from the velocity field.
Definition: transport_dg.cc:671
TransportDG::EqData::set_DG_parameters_boundary
void set_DG_parameters_boundary(Side 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.
Definition: transport_dg.cc:171
index_types.hh
TransportDG::mass_assembly_
GenericAssembly< MassAssemblyDim > * mass_assembly_
general assembly objects, hold assembly objects of appropriate dimension
Definition: transport_dg.hh:421
TransportDG::mass_matrix
std::vector< Mat > mass_matrix
The mass matrix.
Definition: transport_dg.hh:378
TransportDG::calculate_cumulative_balance
void calculate_cumulative_balance()
Definition: transport_dg.cc:650
TransportDG::compute_p0_interpolation
void compute_p0_interpolation()
Compute P0 interpolation of the solution (used in reaction term).
Definition: transport_dg.cc:598
storage.hh
type_record.hh
TransportDG::EqFields::region_id
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: transport_dg.hh:154
ASSERT_PTR_DBG
#define ASSERT_PTR_DBG(ptr)
Definition of assert macro checking non-null pointer (PTR) only for debug mode.
Definition: asserts.hh:340
TransportDG::TransportDG
TransportDG(Mesh &init_mesh, const Input::Record in_rec)
Constructor.
Definition: transport_dg.cc:207
TransportDG::evaluate_time_constraint
bool evaluate_time_constraint(double &)
Definition: transport_dg.hh:277
TransportDG::EqData::ret_vec
std::vector< Vec > ret_vec
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
Definition: transport_dg.hh:214
TransportDG::get_component_vec
Vec get_component_vec(unsigned int sbi)
Return PETSc vector with solution for sbi-th component.
Definition: transport_dg.hh:300
TransportDG::EqData::dg_variant
int dg_variant
DG variant ((non-)symmetric/incomplete.
Definition: transport_dg.hh:205
equation_output.hh
TransportDG::EqData
Definition: transport_dg.hh:161
Distribution
Definition: distribution.hh:50
TransportDG::rhs
std::vector< Vec > rhs
Vector of right hand side.
Definition: transport_dg.hh:372
TransportDG::EqData::ls
LinSys ** ls
Linear algebra system for the transport equation.
Definition: transport_dg.hh:217
TransportDG::output_data
void output_data()
Postprocesses the solution and writes to output file.
Definition: transport_dg.cc:625
Mapping
Definition: transport_dg.hh:68
accessors.hh
TransportDG::eq_fields_
std::shared_ptr< EqFields > eq_fields_
Fields for model parameters.
Definition: transport_dg.hh:359
TransportDG::EqData::balance_
std::shared_ptr< Balance > balance_
Definition: transport_dg.hh:241
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
elements.h
type_generic.hh
TransportDG::eq_fields
Model::ModelEqFields & eq_fields()
Definition: transport_dg.hh:317
DOFHandlerMultiDim
Provides the numbering of the finite element degrees of freedom on the computational mesh.
Definition: dofhandler.hh:151
mpi.h
TransportDG::EqData::dh_
std::shared_ptr< DOFHandlerMultiDim > dh_
Object for distribution of dofs.
Definition: transport_dg.hh:233
TransportDG::initialize
void initialize() override
Definition: transport_dg.cc:247
accessors.hh
EquationOutput
Definition: equation_output.hh:44
TimeGovernor
Basic time management functionality for unsteady (and steady) solvers (class Equation).
Definition: time_governor.hh:317
TransportDG::get_dg_variant_selection_input_type
static const Input::Type::Selection & get_dg_variant_selection_input_type()
Input type for the DG variant selection.
Definition: transport_dg.cc:53
TransportDG::preallocate
void preallocate()
Definition: transport_dg.cc:432
Side
Definition: accessors.hh:390
field_values.hh
OutputTime
The class for outputting data during time.
Definition: output_time.hh:51
AssemblyDG
Definition: transport_dg.hh:58
SourcesAssemblyDG
Definition: assembly_dg.hh:651
TransportDG::allocation_done
bool allocation_done
Indicates whether matrices have been preallocated.
Definition: transport_dg.hh:415
TransportDG::EqFields
Definition: transport_dg.hh:147
TransportDG::EqData::dif_coef
vector< vector< arma::mat33 > > dif_coef
Diffusion coefficients.
Definition: transport_dg.hh:226
TransportDG::balance
std::shared_ptr< Balance > balance() const
Access to balance object of Model.
Definition: transport_dg.hh:313
Input::Type::Selection
Template for classes storing finite set of named values.
Definition: type_selection.hh:65
heat_model.hh
Discontinuous Galerkin method for equation of transport with dispersion.
FiniteElement
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: discrete_space.hh:26
FEValuesBase
Definition: transport_dg.hh:66
BdrConditionAssemblyDG
Definition: assembly_dg.hh:776
TransportDG::registrar
static const int registrar
Registrar of class to factory.
Definition: transport_dg.hh:323
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
TransportDG
Transport with dispersion implemented using discontinuous Galerkin method.
Definition: transport_dg.hh:134
TransportDG::EqFields::fracture_sigma
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
Definition: transport_dg.hh:152
TransportDG::stiffness_matrix
std::vector< Mat > stiffness_matrix
The stiffness matrix.
Definition: transport_dg.hh:375
TransportDG::EqData::dh_p0
std::shared_ptr< DOFHandlerMultiDim > dh_p0
Definition: transport_dg.hh:239
TransportDG::eq_data_
std::shared_ptr< EqData > eq_data_
Data for model parameters.
Definition: transport_dg.hh:362
TransportDG::DGVariant
DGVariant
Definition: transport_dg.hh:244
TransportDG::reg_stat_stream
ofstream reg_stat_stream
Definition: transport_dg.hh:394
Mesh
Definition: mesh.h:355
TransportDG::bdr_cond_assembly_
GenericAssembly< BdrConditionAssemblyDim > * bdr_cond_assembly_
Definition: transport_dg.hh:424
multi_field.hh
TransportDG::EqData::EqData
EqData()
Definition: transport_dg.hh:164
accessors_impl.hh
TransportDG::update_after_reactions
void update_after_reactions(bool solution_changed)
Definition: transport_dg.cc:697
Model
Definition: field_model.hh:338
update_flags.hh
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
TransportDG::EqFields::subdomain
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: transport_dg.hh:155
MultiField
Class for representation of a vector of fields of the same physical quantity.
Definition: multi_field.hh:87
TransportDG::EqData::gamma
std::vector< std::vector< double > > gamma
Penalty parameters.
Definition: transport_dg.hh:202
TransportDG::output_region_statistics
void output_region_statistics()
TransportDG::mass_vec
std::vector< Vec > mass_vec
Mass from previous time instant (necessary when coefficients of mass matrix change in time).
Definition: transport_dg.hh:381
concentration_model.hh
Discontinuous Galerkin method for equation of transport with dispersion.
InitProjectionAssemblyDG
Definition: assembly_dg.hh:998
DiscreteSpace
Definition: discrete_space.hh:41
TransportDG::zero_time_step
void zero_time_step() override
Initialize solution in the zero time.
Definition: transport_dg.cc:399
TransportDG::~TransportDG
~TransportDG() override
Destructor.
Definition: transport_dg.cc:357
TransportDG::EqData::output_vec
std::vector< VectorMPI > output_vec
Vector of solution data.
Definition: transport_dg.hh:236
InitConditionAssemblyDG
Definition: assembly_dg.hh:1097
Field
Class template representing a field with values dependent on: point, element, and region.
Definition: field.hh:93
TransportDG::eq_data
Model::ModelEqData & eq_data()
Definition: transport_dg.hh:319
TransportDG::sources_assembly_
GenericAssembly< SourcesAssemblyDim > * sources_assembly_
Definition: transport_dg.hh:423
TransportDG::init_assembly_
GenericAssemblyBase * init_assembly_
Definition: transport_dg.hh:425
TransportDG::get_input_type
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
Definition: transport_dg.cc:79
Quadrature
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
TransportDG::ret_sources
vector< double > ret_sources
Temporary values of increments due to retardation (e.g. sorption)
Definition: transport_dg.hh:404
mixed.hh
TransportDG::non_symmetric
@ non_symmetric
Definition: transport_dg.hh:246
TransportDG::EqData::elem_anisotropy
double elem_anisotropy(ElementAccessor< 3 > e) const
Compute and return anisotropy of given element.
Definition: transport_dg.cc:155
GenericAssembly
Generic class of assemblation.
Definition: generic_assembly.hh:160
TransportDG::EqFields::output_fields
EquationOutput output_fields
Definition: transport_dg.hh:157
StiffnessAssemblyDG
Definition: assembly_dg.hh:163
field.hh
IntDim
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19