Flow123d  release_3.0.0-880-gc768b7a
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 #include <boost/exception/info.hpp> // for operator<<, error_info...
26 #include <string> // for operator<<
27 #include <vector> // for vector
28 #include <armadillo>
29 #include "fem/update_flags.hh" // for operator|
30 #include "fem/mapping_p1.hh"
31 #include "fields/field_values.hh" // for FieldValue<>::Scalar
32 #include "fields/field.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 "mesh/long_idx.hh" // for LongIdx
44 #include "mesh/accessors.hh" // for ElementAccessor
45 #include "mesh/elements.h" // for Element::dim, Element:...
46 #include "mesh/neighbours.h" // for Neighbour::element
47 #include "mesh/side_impl.hh" // for Side::cond, Side::cond...
48 #include "mesh/sides.h" // for SideIter
49 #include "mpi.h" // for MPI_Comm_rank
50 #include "petscmat.h" // for Mat, MatDestroy
51 #include "petscvec.h" // for Vec, VecDestroy, VecSc...
52 #include "transport/concentration_model.hh" // for ConcentrationTransport...
53 #include "transport/heat_model.hh" // for HeatTransferModel, Hea...
54 
55 class DiscreteSpace;
56 class Distribution;
57 class OutputTime;
58 class DOFHandlerMultiDim;
59 template<unsigned int dim, unsigned int spacedim> class FEValuesBase;
60 template<unsigned int dim> class FiniteElement;
61 template<unsigned int dim, unsigned int spacedim> class Mapping;
62 template<unsigned int dim> class Quadrature;
63 namespace Input { namespace Type { class Selection; } }
64 
65 
66 
67 /**
68  * Auxiliary container class for Finite element and related objects of all dimensions.
69  * Its purpose is to provide templated access to these objects, applicable in
70  * the assembling methods.
71  */
72 class FEObjects {
73 public:
74 
75  FEObjects(Mesh *mesh_, unsigned int fe_order);
76  ~FEObjects();
77 
78  template<unsigned int dim>
79  inline FiniteElement<dim> *fe();
80 
81  template<unsigned int dim>
82  inline FiniteElement<dim> *fe_rt();
83 
84  template<unsigned int dim>
85  inline Quadrature<dim> *q();
86 
87  template<unsigned int dim>
88  inline MappingP1<dim,3> *mapping();
89 
90  inline std::shared_ptr<DOFHandlerMultiDim> dh();
91 
92 private:
93 
94  /// Finite elements for the solution of the advection-diffusion equation.
99 
100  /// Finite elements for the water velocity field.
104 
105  /// Quadratures used in assembling methods.
110 
111  /// Auxiliary mappings of reference elements.
115 
116  std::shared_ptr<DiscreteSpace> ds_;
117 
118  /// Object for distribution of dofs.
119  std::shared_ptr<DOFHandlerMultiDim> dh_;
120 };
121 
122 
123 
124 /**
125  * @brief Transport with dispersion implemented using discontinuous Galerkin method.
126  *
127  * TransportDG implements the discontinuous Galerkin method for the transport and diffusion of substances.
128  * The concentration @f$ c_i ~[kg/m^3]@f$ of the i-th substance is governed by the advection-diffusion equation
129  * @f[
130  * \partial_t c_i + \mathbf v\cdot\nabla c_i - \mathrm{div}(D\nabla c_i) = F \mbox{ in }\Omega^d,
131  * @f]
132  * where @f$\mathbf v@f$ is the fluid velocity and @f$\Omega^d@f$ the @f$d@f$-dimensional domain, respectively.
133  * The hydrodynamic dispersivity tensor @f$\mathbf D ~[m^2/s]@f$ is given by:
134  * @f[
135  * \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).
136  * @f]
137  * 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.
138  *
139  * 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$.
140  *
141  * 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$.
142  * We prescribe the following boundary conditions:
143  * @f{eqnarray*}{
144  * c_i^d &= c_{iD}^d &\mbox{ on }\Gamma^d_D \mbox{ (Dirichlet)},\\
145  * \mathbf D^d\nabla c_i^d\cdot\mathbf n &= 0 &\mbox{ on }\Gamma^d_N \mbox{ (Neumann)},
146  * @f}
147  * The transfer of mass through fractures is described by the transmission conditions on @f$\Gamma^d_F@f$:
148  * @f[
149  * -\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
150  * F^{d-1} = (\sigma + |\mathbf v^d\cdot\mathbf n|)(c_i^d-c_i^{d-1}).
151  * @f]
152  * Here @f$\mathbf n@f$ stands for the unit outward normal vector to @f$\partial\Omega^d@f$.
153  * The coefficient @f$\sigma@f$ determines the transfer of mass through fractures due to diffusion.
154  *
155  * @ingroup transport_mod
156  *
157  */
158 template<class Model>
159 class TransportDG : public Model
160 {
161 public:
162 
163  class EqData : public Model::ModelEqData {
164  public:
165 
166  EqData();
167 
168  MultiField<3, FieldValue<3>::Scalar> fracture_sigma; ///< Transition parameter for diffusive transfer on fractures (for each substance).
169  MultiField<3, FieldValue<3>::Scalar> dg_penalty; ///< Penalty enforcing inter-element continuity of solution (for each substance).
172 
174 
175  };
176 
177 
178 
179  enum DGVariant {
180  // Non-symmetric weighted interior penalty DG
181  non_symmetric = -1,
182 
183  // Incomplete weighted interior penalty DG
184  incomplete = 0,
185 
186  // Symmetric weighted interior penalty DG
187  symmetric = 1
188  };
189 
190  /**
191  * @brief Constructor.
192  * @param init_mesh computational mesh
193  * @param in_rec input record
194  */
195  TransportDG(Mesh &init_mesh, const Input::Record in_rec);
196  /**
197 
198  * @brief Declare input record type for the equation TransportDG.
199  */
200  static const Input::Type::Record & get_input_type();
201 
202  /**
203  * @brief Input type for the DG variant selection.
204  */
205  static const Input::Type::Selection & get_dg_variant_selection_input_type();
206 
207  /**
208  * @brief Initialize solution in the zero time.
209  */
210  void zero_time_step() override;
211 
212  bool evaluate_time_constraint(double &time_constraint)
213  { return false; }
214 
215  /**
216  * @brief Computes the solution in one time instant.
217  */
218  void update_solution() override;
219 
220  /**
221  * @brief Postprocesses the solution and writes to output file.
222  */
223  void output_data();
224 
225  /**
226  * @brief Destructor.
227  */
228  ~TransportDG();
229 
230  void initialize() override;
231 
232  void calculate_cumulative_balance();
233 
234  const Vec &get_solution(unsigned int sbi)
235  { return ls[sbi]->get_solution(); }
236 
238  { return solution_elem_; }
239 
240  void calculate_concentration_matrix();
241 
242  void update_after_reactions(bool solution_changed);
243 
244  void get_par_info(LongIdx * &el_4_loc, Distribution * &el_ds);
245 
246  LongIdx *get_row_4_el();
247 
248 
249 
250 
251 private:
252  /// Registrar of class to factory
253  static const int registrar;
254 
255  inline typename Model::ModelEqData &data() { return data_; }
256 
257  void preallocate();
258 
259  /**
260  * @brief Assembles the mass matrix.
261  *
262  * The routine just calls templated method assemble_mass_matrix() for each
263  * space dimension.
264  */
265  void assemble_mass_matrix();
266 
267  /**
268  * @brief Assembles the mass matrix for the given dimension.
269  */
270  template<unsigned int dim>
271  void assemble_mass_matrix();
272 
273  /**
274  * @brief Assembles the stiffness matrix.
275  *
276  * This routine just calls assemble_volume_integrals(), assemble_fluxes_boundary(),
277  * assemble_fluxes_element_element() and assemble_fluxes_element_side() for each
278  * space dimension.
279  */
280  void assemble_stiffness_matrix();
281 
282  /**
283  * @brief Assembles the volume integrals into the stiffness matrix.
284  */
285  template<unsigned int dim>
286  void assemble_volume_integrals();
287 
288  /**
289  * @brief Assembles the right hand side due to volume sources.
290  *
291  * This method just calls set_sources() for each space dimension.
292  */
293  void set_sources();
294 
295  /**
296  * @brief Assembles the right hand side vector due to volume sources.
297  */
298  template<unsigned int dim>
299  void set_sources();
300 
301  /**
302  * @brief Assembles the fluxes on the boundary.
303  */
304  template<unsigned int dim>
305  void assemble_fluxes_boundary();
306 
307  /**
308  * @brief Assembles the fluxes between elements of the same dimension.
309  */
310  template<unsigned int dim>
311  void assemble_fluxes_element_element();
312 
313  /**
314  * @brief Assembles the fluxes between elements of different dimensions.
315  */
316  template<unsigned int dim>
317  void assemble_fluxes_element_side();
318 
319 
320  /**
321  * @brief Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
322  *
323  * The routine just calls templated method set_boundary_condition() for each space dimension.
324  */
325  void set_boundary_conditions();
326 
327  /**
328  * @brief Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions
329  * for a given space dimension.
330  */
331  template<unsigned int dim>
332  void set_boundary_conditions();
333 
334  /**
335  * @brief Calculates the velocity field on a given @p dim dimensional cell.
336  *
337  * @param cell The cell.
338  * @param velocity The computed velocity field (at quadrature points).
339  * @param fv The FEValues class providing the quadrature points
340  * and the shape functions for velocity.
341  */
342  template<unsigned int dim>
343  void calculate_velocity(const ElementAccessor<3> &cell, std::vector<arma::vec3> &velocity, FEValuesBase<dim,3> &fv);
344 
345  /**
346  * @brief Calculates the dispersivity (diffusivity) tensor from the velocity field.
347  *
348  * @param K The computed dispersivity tensor.
349  * @param velocity The velocity field (at quadrature points).
350  * @param Dm Molecular diffusivities.
351  * @param alphaL Longitudal dispersivities.
352  * @param alphaT Transversal dispersivities.
353  * @param porosity Porosities.
354  * @param cross_cut Cross-cuts of higher dimension.
355  */
356 // void calculate_dispersivity_tensor(arma::mat33 &K, const arma::vec3 &velocity,
357 // double Dm, double alphaL, double alphaT, double porosity,
358 // double cross_cut);
359 
360  /**
361  * @brief Sets up some parameters of the DG method for two sides of an edge.
362  *
363  * @param edg The edge.
364  * @param s1 Side 1.
365  * @param s2 Side 2.
366  * @param K_size Size of vector of tensors K.
367  * @param K1 Dispersivity tensors on side s1 (in quadrature points).
368  * @param K2 Dispersivity tensors on side s2 (in quadrature points).
369  * @param normal_vector Normal vector to side 0 of the neighbour
370  * (assumed constant along the side).
371  * @param alpha1, alpha2 Penalty parameter that influences the continuity
372  * of the solution (large value=more continuity).
373  * @param gamma Computed penalty parameters.
374  * @param omega Computed weights.
375  * @param transport_flux Computed flux from side s1 to side s2.
376  */
377  void set_DG_parameters_edge(const Edge &edg,
378  const int s1,
379  const int s2,
380  const int K_size,
381  const std::vector<arma::mat33> &K1,
382  const std::vector<arma::mat33> &K2,
383  const std::vector<double> &fluxes,
384  const arma::vec3 &normal_vector,
385  const double alpha1,
386  const double alpha2,
387  double &gamma,
388  double *omega,
389  double &transport_flux);
390 
391  /**
392  * @brief Sets up parameters of the DG method on a given boundary edge.
393  *
394  * Assumption is that the edge consists of only 1 side.
395  * @param side The boundary side.
396  * @param K_size Size of vector of tensors K.
397  * @param K Dispersivity tensor.
398  * @param ad_vector Advection vector.
399  * @param normal_vector Normal vector (assumed constant along the edge).
400  * @param alpha Penalty parameter that influences the continuity
401  * of the solution (large value=more continuity).
402  * @param gamma Computed penalty parameters.
403  */
404  void set_DG_parameters_boundary(const SideIter side,
405  const int K_size,
406  const std::vector<arma::mat33> &K,
407  const double flux,
408  const arma::vec3 &normal_vector,
409  const double alpha,
410  double &gamma);
411 
412 
413  /**
414  * @brief Sets the initial condition.
415  */
416  void set_initial_condition();
417 
418  /**
419  * @brief Assembles the auxiliary linear system to calculate the initial solution
420  * as L^2-projection of the prescribed initial condition.
421  */
422  template<unsigned int dim>
423  void prepare_initial_condition();
424 
425 
426 
427  /// @name Physical parameters
428  // @{
429 
430  /// Field data for model parameters.
432 
433  // @}
434 
435 
436  /// @name Parameters of the numerical method
437  // @{
438 
439  /// Finite element objects
441 
442  /// Penalty parameters.
444 
445  /// DG variant ((non-)symmetric/incomplete
447 
448  /// Polynomial order of finite elements.
449  unsigned int dg_order;
450 
451  // @}
452 
453 
454 
455  /// @name Solution of algebraic system
456  // @{
457 
458  /// Vector of right hand side.
460 
461  /// The stiffness matrix.
463 
464  /// The mass matrix.
466 
467  /// Mass from previous time instant (necessary when coefficients of mass matrix change in time).
469 
470  /// Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
472 
473  /// Linear algebra system for the transport equation.
475 
476  /// Linear algebra system for the time derivative (actually it is used only for handling the matrix structures).
478 
479  /// Element averages of solution (the array is passed to reactions in operator splitting).
480  double **solution_elem_;
481 
482  // @}
483 
484 
485  /// @name Output to file
486  // @{
487 
488  /// Array for storing the output solution data.
489  //vector<double*> output_solution;
490 
491  /// Vector of solution data.
493 
494  /// Record with input specification.
496 
497 
498  // @}
499 
500 
501  /// @name Auxiliary fields used during assembly
502  // @{
503 
504  /// Mass matrix coefficients.
506  /// Retardation coefficient due to sorption.
508  /// Temporary values of increments due to retardation (e.g. sorption)
510  /// Advection coefficients.
512  /// Diffusion coefficients.
514  /// Advection coefficients on edges.
516  /// Diffusion coefficients on edges.
518 
519  // @}
520 
521 
522 
523 
524  /// @name Other
525  // @{
526 
527  /// Indicates whether matrices have been preallocated.
529 
530  // @}
531 };
532 
533 
534 
535 
536 
537 
538 #endif /* TRANSPORT_DG_HH_ */
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
Input::Record input_rec
Record with input specification.
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
FiniteElement< 2 > * fe_rt2_
vector< double > mm_coef
Mass matrix coefficients.
MappingP1< 3, 3 > * map3_
Transport with dispersion implemented using discontinuous Galerkin method.
Quadrature< 1 > * q1_
Abstract linear system class.
Definition: balance.hh:35
Wrappers for linear systems based on MPIAIJ and MATIS format.
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:83
vector< VectorMPI > output_vec
Array for storing the output solution data.
int dg_variant
DG variant ((non-)symmetric/incomplete.
Definition: mesh.h:80
vector< vector< arma::vec3 > > ad_coef
Advection coefficients.
Definition: edges.h:26
vector< vector< vector< arma::mat33 > > > dif_coef_edg
Diffusion coefficients on edges.
Discontinuous Galerkin method for equation of transport with dispersion.
const Vec & get_solution(unsigned int sbi)
vector< double > ret_sources_prev
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
Quadrature< 2 > * q2_
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
EquationOutput output_fields
std::vector< Mat > mass_matrix
The mass matrix.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:35
vector< vector< double > > ret_coef
Retardation coefficient due to sorption.
vector< vector< vector< arma::vec3 > > > ad_coef_edg
Advection coefficients on edges.
double ** get_concentration_matrix()
FiniteElement< 0 > * fe0_
Finite elements for the solution of the advection-diffusion equation.
Definition: transport_dg.hh:95
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< 2 > * fe2_
Definition: transport_dg.hh:97
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:38
Field< 3, FieldValue< 3 >::Scalar > subdomain
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
Quadrature< 3 > * q3_
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:152
std::shared_ptr< DiscreteSpace > ds_
std::vector< Vec > rhs
Vector of right hand side.
bool evaluate_time_constraint(double &time_constraint)
Field< 3, FieldValue< 3 >::Scalar > region_id
unsigned int dg_order
Polynomial order of finite elements.
The class for outputting data during time.
Definition: output_time.hh:50
FiniteElement< 1 > * fe_rt1_
Finite elements for the water velocity field.
MappingP1< 2, 3 > * map2_
double ** solution_elem_
Element averages of solution (the array is passed to reactions in operator splitting).
Model::ModelEqData & data()
std::vector< Mat > stiffness_matrix
The stiffness matrix.
EqData data_
Field data for model parameters.
bool allocation_done
Indicates whether matrices have been preallocated.
std::vector< std::vector< double > > gamma
Penalty parameters.
FiniteElement< 1 > * fe1_
Definition: transport_dg.hh:96
std::vector< Vec > ret_vec
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
FiniteElement< 3 > * fe3_
Definition: transport_dg.hh:98
FiniteElement< 3 > * fe_rt3_
Discontinuous Galerkin method for equation of transport with dispersion.
std::vector< Vec > mass_vec
Mass from previous time instant (necessary when coefficients of mass matrix change in time)...
Record type proxy class.
Definition: type_record.hh:182
Class for representation of a vector of fields of the same physical quantity.
Definition: multi_field.hh:87
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Base class for FEValues and FESideValues.
Definition: fe_values.hh:37
MappingP1< 1, 3 > * map1_
Auxiliary mappings of reference elements.
LinSys ** ls
Linear algebra system for the transport equation.
LinSys ** ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
Template for classes storing finite set of named values.
Quadrature< 0 > * q0_
Quadratures used in assembling methods.
std::shared_ptr< DOFHandlerMultiDim > dh_
Object for distribution of dofs.