Flow123d  DF_patch_fe_data_tables-65e7202
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> class FiniteElement;
67 template<unsigned int dim, unsigned int spacedim> class Mapping;
68 class Quadrature;
69 namespace Input { namespace Type { class Selection; } }
70 class ElementCacheMap;
71 
72 /*class FEObjects {
73 public:
74 
75  inline FEObjects(Mesh *mesh_, unsigned int fe_order)
76  {
77  fe = MixedPtr<FE_P_disc>(fe_order);
78  fe_rt = MixedPtr<FE_RT0>();
79  q = MixedPtr<QGauss>(2*fe_order);
80 
81  auto ds = std::make_shared<EqualOrderDiscreteSpace>(mesh_, fe);
82  dh = std::make_shared<DOFHandlerMultiDim>(*mesh_);
83  dh->distribute_dofs(ds_);
84  }
85 
86  ~FEObjects();
87 
88  MixedPtr<FiniteElement> fe;
89  MixedPtr<FiniteElement> fe_rt;
90  MixedPtr<Quadrature> q;
91 
92  /// Object for distribution of dofs.
93  std::shared_ptr<DOFHandlerMultiDim> dh;
94 };*/
95 
96 
97 
98 /**
99  * @brief Transport with dispersion implemented using discontinuous Galerkin method.
100  *
101  * TransportDG implements the discontinuous Galerkin method for the transport and diffusion of substances.
102  * The concentration @f$ c_i ~[kg/m^3]@f$ of the i-th substance is governed by the advection-diffusion equation
103  * @f[
104  * \partial_t c_i + \mathbf v\cdot\nabla c_i - \mathrm{div}(D\nabla c_i) = F \mbox{ in }\Omega^d,
105  * @f]
106  * where @f$\mathbf v@f$ is the fluid velocity and @f$\Omega^d@f$ the @f$d@f$-dimensional domain, respectively.
107  * The hydrodynamic dispersivity tensor @f$\mathbf D ~[m^2/s]@f$ is given by:
108  * @f[
109  * \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).
110  * @f]
111  * 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.
112  *
113  * 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$.
114  *
115  * 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$.
116  * We prescribe the following boundary conditions:
117  * @f{eqnarray*}{
118  * c_i^d &= c_{iD}^d &\mbox{ on }\Gamma^d_D \mbox{ (Dirichlet)},\\
119  * \mathbf D^d\nabla c_i^d\cdot\mathbf n &= 0 &\mbox{ on }\Gamma^d_N \mbox{ (Neumann)},
120  * @f}
121  * The transfer of mass through fractures is described by the transmission conditions on @f$\Gamma^d_F@f$:
122  * @f[
123  * -\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
124  * F^{d-1} = (\sigma + |\mathbf v^d\cdot\mathbf n|)(c_i^d-c_i^{d-1}).
125  * @f]
126  * Here @f$\mathbf n@f$ stands for the unit outward normal vector to @f$\partial\Omega^d@f$.
127  * The coefficient @f$\sigma@f$ determines the transfer of mass through fractures due to diffusion.
128  *
129  * @ingroup transport_mod
130  *
131  */
132 template<class Model>
133 class TransportDG : public Model
134 {
135 public:
136 
137  template<unsigned int dim> using MassAssemblyDim = MassAssemblyDG<dim, Model>;
138  template<unsigned int dim> using StiffnessAssemblyDim = StiffnessAssemblyDG<dim, Model>;
139  template<unsigned int dim> using SourcesAssemblyDim = SourcesAssemblyDG<dim, Model>;
143 
145 
146  class EqFields : public Model::ModelEqFields {
147  public:
148 
149  EqFields();
150 
151  MultiField<3, FieldValue<3>::Scalar> fracture_sigma; ///< Transition parameter for diffusive transfer on fractures (for each substance).
152  MultiField<3, FieldValue<3>::Scalar> dg_penalty; ///< Penalty enforcing inter-element continuity of solution (for each substance).
155 
157  };
158 
159 
160  class EqData : public Model::ModelEqData {
161  public:
162 
163  EqData() : ls(nullptr) {}
164 
165  inline void set_time_governor(TimeGovernor *time) {
166  ASSERT_PTR(time);
167  this->time_ = time;
168  }
169 
170  /// Getter of dg_order
171  inline unsigned int quad_order() const {
172  return dg_order;
173  }
174 
175  /// @name Parameters of the numerical method
176  // @{
177 
178  /// DG variant ((non-)symmetric/incomplete
180 
181  /// Polynomial order of finite elements.
182  unsigned int dg_order;
183 
184  // @}
185 
186 
187  /// Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
189 
190  /// Linear algebra system for the transport equation.
192 
193  /// Linear algebra system for the time derivative (actually it is used only for handling the matrix structures).
195 
196  /// @name Auxiliary fields used during assembly
197  // @{
198 
199  /// Maximal number of edge sides (evaluate from dim 1,2,3)
200  unsigned int max_edg_sides;
201 
202  // @}
203 
204  /// Object for distribution of dofs.
205  std::shared_ptr<DOFHandlerMultiDim> dh_;
206 
207  /// Vector of solution data.
209 
211  std::shared_ptr<DOFHandlerMultiDim> dh_p0;
213  std::shared_ptr<Balance> balance_;
214  };
215 
216  enum DGVariant {
217  // Non-symmetric weighted interior penalty DG
219 
220  // Incomplete weighted interior penalty DG
222 
223  // Symmetric weighted interior penalty DG
224  symmetric = 1
225  };
226 
227  /**
228  * @brief Constructor.
229  * @param init_mesh computational mesh
230  * @param in_rec input record
231  */
232  TransportDG(Mesh &init_mesh, const Input::Record in_rec);
233  /**
234 
235  * @brief Declare input record type for the equation TransportDG.
236  */
237  static const Input::Type::Record & get_input_type();
238 
239  /**
240  * @brief Input type for the DG variant selection.
241  */
243 
244  /**
245  * @brief Initialize solution in the zero time.
246  */
247  void zero_time_step() override;
248 
250  { return false; }
251 
252  /**
253  * @brief Computes the solution in one time instant.
254  */
255  void update_solution() override;
256 
257  /**
258  * @brief Postprocesses the solution and writes to output file.
259  */
260  void output_data();
261 
262  /**
263  * @brief Destructor.
264  */
265  ~TransportDG() override;
266 
267  void initialize() override;
268 
270 
271  /// Return PETSc vector with solution for sbi-th component.
272  Vec get_component_vec(unsigned int sbi)
273  { return eq_data_->ls[sbi]->get_solution(); }
274 
275  /// Getter for P0 interpolation by FieldFE.
277  { return eq_data_->conc_fe;}
278 
279  /// Compute P0 interpolation of the solution (used in reaction term).
281 
282  void update_after_reactions(bool solution_changed);
283 
284  /// Access to balance object of Model
285  inline std::shared_ptr<Balance> balance() const {
286  return Model::balance_;
287  }
288 
289  inline typename Model::ModelEqFields &eq_fields() { return *eq_fields_; }
290 
291  inline typename Model::ModelEqData &eq_data() { return *eq_data_; }
292 
293 private:
294  /// Registrar of class to factory
295  static const int registrar;
296 
297  void preallocate();
298 
299  /**
300  * @brief Calculates the dispersivity (diffusivity) tensor from the velocity field.
301  *
302  * @param K The computed dispersivity tensor.
303  * @param velocity The velocity field (at quadrature points).
304  * @param Dm Molecular diffusivities.
305  * @param alphaL Longitudal dispersivities.
306  * @param alphaT Transversal dispersivities.
307  * @param porosity Porosities.
308  * @param cross_cut Cross-cuts of higher dimension.
309  */
310 // void calculate_dispersivity_tensor(arma::mat33 &K, const arma::vec3 &velocity,
311 // double Dm, double alphaL, double alphaT, double porosity,
312 // double cross_cut);
313 
314  /**
315  * @brief Sets the initial condition.
316  *
317  * This method just calls AssemblyDG::prepare_initial_condition() for each elements.
318  */
319  void set_initial_condition();
320 
321 
322 
324 
325 
326 
327  /// @name Physical parameters
328  // @{
329 
330  /// Fields for model parameters.
331  std::shared_ptr<EqFields> eq_fields_;
332 
333  /// Data for model parameters.
334  std::shared_ptr<EqData> eq_data_;
335 
336  // @}
337 
338 
339 
340  /// @name Solution of algebraic system
341  // @{
342 
343  /// Vector of right hand side.
345 
346  /// The stiffness matrix.
348 
349  /// The mass matrix.
351 
352  /// Mass from previous time instant (necessary when coefficients of mass matrix change in time).
354  // @}
355 
356 
357  /// @name Output to file
358  // @{
359 
360  /// Array for storing the output solution data.
361  //vector<double*> output_solution;
362 
363  /// Record with input specification.
365 
366  ofstream reg_stat_stream;
367 
368 
369  // @}
370 
371 
372  /// @name Auxiliary fields used during assembly
373  // @{
374 
375  /// Temporary values of increments due to retardation (e.g. sorption)
377 
378  // @}
379 
380 
381 
382 
383  /// @name Other
384  // @{
385 
386  /// Indicates whether matrices have been preallocated.
388 
390  // @}
391 
392  /// general assembly objects, hold assembly objects of appropriate dimension
398 };
399 
400 
401 
402 
403 
404 
405 #endif /* TRANSPORT_DG_HH_ */
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR) only for debug mode.
Definition: asserts.hh:341
Provides the numbering of the finite element degrees of freedom on the computational mesh.
Definition: dofhandler.hh:151
Directing class of FieldValueCache.
Class template representing a field with values dependent on: point, element, and region.
Definition: field.hh:92
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Generic class of assemblation.
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
Record type proxy class.
Definition: type_record.hh:182
Template for classes storing finite set of named values.
Definition: mesh.h:362
Class for representation of a vector of fields of the same physical quantity.
Definition: multi_field.hh:87
The class for outputting data during time.
Definition: output_time.hh:51
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Basic time management functionality for unsteady (and steady) solvers (class Equation).
FieldFEScalarVec conc_fe
std::shared_ptr< DOFHandlerMultiDim > dh_
Object for distribution of dofs.
std::vector< VectorMPI > output_vec
Vector of solution data.
LinSys ** ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
std::shared_ptr< DOFHandlerMultiDim > dh_p0
int dg_variant
DG variant ((non-)symmetric/incomplete.
TimeGovernor * time_
void set_time_governor(TimeGovernor *time)
std::shared_ptr< Balance > balance_
unsigned int dg_order
Polynomial order of finite elements.
LinSys ** ls
Linear algebra system for the transport equation.
std::vector< Vec > ret_vec
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
unsigned int quad_order() const
Getter of dg_order.
unsigned int max_edg_sides
Maximal number of edge sides (evaluate from dim 1,2,3)
Field< 3, FieldValue< 3 >::Scalar > region_id
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
EquationOutput output_fields
Field< 3, FieldValue< 3 >::Scalar > subdomain
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
Transport with dispersion implemented using discontinuous Galerkin method.
std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > FieldFEScalarVec
void update_solution() override
Computes the solution in one time instant.
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
Definition: transport_dg.cc:78
std::vector< Vec > mass_vec
Mass from previous time instant (necessary when coefficients of mass matrix change in time).
bool allocation_done
Indicates whether matrices have been preallocated.
void output_region_statistics()
Vec get_component_vec(unsigned int sbi)
Return PETSc vector with solution for sbi-th component.
Input::Record input_rec
Array for storing the output solution data.
void zero_time_step() override
Initialize solution in the zero time.
std::shared_ptr< EqFields > eq_fields_
Fields for model parameters.
void initialize() override
GenericAssembly< MassAssemblyDim > * mass_assembly_
general assembly objects, hold assembly objects of appropriate dimension
static const Input::Type::Selection & get_dg_variant_selection_input_type()
Input type for the DG variant selection.
Definition: transport_dg.cc:52
GenericAssemblyBase * init_assembly_
vector< double > ret_sources
Temporary values of increments due to retardation (e.g. sorption)
FieldFEScalarVec & get_p0_interpolation()
Getter for P0 interpolation by FieldFE.
GenericAssembly< SourcesAssemblyDim > * sources_assembly_
GenericAssembly< BdrConditionAssemblyDim > * bdr_cond_assembly_
Model::ModelEqData & eq_data()
TransportDG(Mesh &init_mesh, const Input::Record in_rec)
Constructor.
~TransportDG() override
Destructor.
std::vector< Vec > rhs
Vector of right hand side.
bool init_projection
GenericAssembly< StiffnessAssemblyDim > * stiffness_assembly_
vector< double > ret_sources_prev
std::vector< Mat > stiffness_matrix
The stiffness matrix.
void compute_p0_interpolation()
Compute P0 interpolation of the solution (used in reaction term).
void output_data()
Postprocesses the solution and writes to output file.
static const int registrar
Registrar of class to factory.
ofstream reg_stat_stream
void calculate_cumulative_balance()
std::shared_ptr< Balance > balance() const
Access to balance object of Model.
std::vector< Mat > mass_matrix
The mass matrix.
void preallocate()
void set_initial_condition()
Calculates the dispersivity (diffusivity) tensor from the velocity field.
std::shared_ptr< EqData > eq_data_
Data for model parameters.
void update_after_reactions(bool solution_changed)
Model::ModelEqFields & eq_fields()
bool evaluate_time_constraint(double &)
Discontinuous Galerkin method for equation of transport with dispersion.
Discontinuous Galerkin method for equation of transport with dispersion.
Wrappers for linear systems based on MPIAIJ and MATIS format.
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19
Abstract linear system class.
Definition: balance.hh:40
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.