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