Flow123d  jenkins-Flow123d-linux-release-multijob-198
transport_dg.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Discontinuous Galerkin method for equation of transport with dispersion.
27  * @author Jan Stebel
28  */
29 
30 #ifndef TRANSPORT_DG_HH_
31 #define TRANSPORT_DG_HH_
32 
34 #include "la/linsys.hh"
35 #include "flow/mh_dofhandler.hh"
36 
37 class Distribution;
38 class OutputTime;
39 class DOFHandlerMultiDim;
40 template<unsigned int dim, unsigned int spacedim> class FEValuesBase;
41 template<unsigned int dim, unsigned int spacedim> class FiniteElement;
42 template<unsigned int dim, unsigned int spacedim> class Mapping;
43 template<unsigned int dim> class Quadrature;
44 
45 
46 
47 /**
48  * Auxiliary container class for Finite element and related objects of all dimensions.
49  * Its purpose is to provide templated access to these objects, applicable in
50  * the assembling methods.
51  */
52 class FEObjects {
53 public:
54 
55  FEObjects(Mesh *mesh_, unsigned int fe_order);
56  ~FEObjects();
57 
58  template<unsigned int dim>
59  inline FiniteElement<dim,3> *fe();
60 
61  template<unsigned int dim>
62  inline FiniteElement<dim,3> *fe_rt();
63 
64  template<unsigned int dim>
65  inline Quadrature<dim> *q();
66 
67  template<unsigned int dim>
68  inline Mapping<dim,3> *mapping();
69 
70  inline DOFHandlerMultiDim *dh();
71 
72 private:
73 
74  /// Finite elements for the solution of the advection-diffusion equation.
78 
79  /// Finite elements for the water velocity field.
83 
84  /// Quadratures used in assembling methods.
89 
90  /// Auxiliary mappings of reference elements.
95 
96  /// Object for distribution of dofs.
98 };
99 
100 
101 
102 /**
103  * @brief Transport with dispersion implemented using discontinuous Galerkin method.
104  *
105  * TransportDG implements the discontinuous Galerkin method for the transport and diffusion of substances.
106  * The concentration @f$ c_i ~[kg/m^3]@f$ of the i-th substance is governed by the advection-diffusion equation
107  * @f[
108  * \partial_t c_i + \mathbf v\cdot\nabla c_i - \mathrm{div}(D\nabla c_i) = F \mbox{ in }\Omega^d,
109  * @f]
110  * where @f$\mathbf v@f$ is the fluid velocity and @f$\Omega^d@f$ the @f$d@f$-dimensional domain, respectively.
111  * The hydrodynamic dispersivity tensor @f$\mathbf D ~[m^2/s]@f$ is given by:
112  * @f[
113  * \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).
114  * @f]
115  * 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.
116  *
117  * 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$.
118  *
119  * 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$.
120  * We prescribe the following boundary conditions:
121  * @f{eqnarray*}{
122  * c_i^d &= c_{iD}^d &\mbox{ on }\Gamma^d_D \mbox{ (Dirichlet)},\\
123  * \mathbf D^d\nabla c_i^d\cdot\mathbf n &= 0 &\mbox{ on }\Gamma^d_N \mbox{ (Neumann)},
124  * @f}
125  * The transfer of mass through fractures is described by the transmission conditions on @f$\Gamma^d_F@f$:
126  * @f[
127  * -\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
128  * F^{d-1} = (\sigma + |\mathbf v^d\cdot\mathbf n|)(c_i^d-c_i^{d-1}).
129  * @f]
130  * Here @f$\mathbf n@f$ stands for the unit outward normal vector to @f$\partial\Omega^d@f$.
131  * The coefficient @f$\sigma@f$ determines the transfer of mass through fractures due to diffusion.
132  *
133  * @ingroup transport_mod
134  *
135  */
136 template<class Model>
137 class TransportDG : public TransportBase, public Model
138 {
139 public:
140 
141  class EqData : public Model::ModelEqData {
142  public:
143 
144  enum BC_Type {
150  };
152 
154 
155  EqData();
156 
157  Field<3, FieldValue<3>::Vector> fracture_sigma; ///< Transition parameter for diffusive transfer on fractures (for each substance).
158  Field<3, FieldValue<3>::Vector> dg_penalty; ///< Penalty enforcing inter-element continuity of solution (for each substance).
159 
160  BCField<3, FieldValue<3>::EnumVector > bc_type; ///< Type of boundary condition (see also BC_Type)
161  BCField<3, FieldValue<3>::Vector > bc_flux; ///< Flux in Neumann or Robin b.c.
162  BCField<3, FieldValue<3>::Vector > bc_robin_sigma; ///< Transition coefficient in Robin b.c.
163 
165 
166  };
167 
168 
169 
170  enum DGVariant {
171  // Non-symmetric weighted interior penalty DG
173 
174  // Incomplete weighted interior penalty DG
176 
177  // Symmetric weighted interior penalty DG
179  };
180 
181  /**
182  * @brief Constructor.
183  * @param init_mesh computational mesh
184  * @param in_rec input record
185  */
186  TransportDG(Mesh &init_mesh, const Input::Record &in_rec);
187  /**
188 
189  * @brief Declare input record type for the equation TransportDG.
190  */
192 
193  /**
194  * @brief Input type for the DG variant selection.
195  */
197 
198  /**
199  * @brief Initialize solution in the zero time.
200  */
201  void zero_time_step() override;
202 
203  /**
204  * @brief Computes the solution in one time instant.
205  */
206  void update_solution() override;
207 
208  /**
209  * @brief Updates the velocity field which determines some coefficients of the transport equation.
210  *
211  * @param dh mixed hybrid dof handler
212  *
213  * (So far it does not work since the flow module returns a vector of zeros.)
214  * @param velocity_vector Input array of velocity values.
215  */
216  virtual void set_velocity_field(const MH_DofHandler &dh);
217 
218  /**
219  * @brief Postprocesses the solution and writes to output file.
220  */
221  void output_data();
222 
223  /**
224  * @brief Getter for field data.
225  */
226  virtual EqData *get_data() { return &data_; }
227 
228  /**
229  * @brief Destructor.
230  */
231  ~TransportDG();
232 
233 private:
234 
235  inline typename Model::ModelEqData &data() { return data_; }
236 
237  void output_vector_gather();
238 
239  void preallocate();
240 
241  /**
242  * @brief Assembles the mass matrix.
243  *
244  * The routine just calls templated method assemble_mass_matrix() for each
245  * space dimension.
246  */
247  void assemble_mass_matrix();
248 
249  /**
250  * @brief Assembles the mass matrix for the given dimension.
251  */
252  template<unsigned int dim>
253  void assemble_mass_matrix();
254 
255  /**
256  * @brief Assembles the stiffness matrix.
257  *
258  * This routine just calls assemble_volume_integrals(), assemble_fluxes_boundary(),
259  * assemble_fluxes_element_element() and assemble_fluxes_element_side() for each
260  * space dimension.
261  */
263 
264  /**
265  * @brief Assembles the volume integrals into the stiffness matrix.
266  */
267  template<unsigned int dim>
269 
270  /**
271  * @brief Assembles the right hand side due to volume sources.
272  *
273  * This method just calls set_sources() for each space dimension.
274  */
275  void set_sources();
276 
277  /**
278  * @brief Assembles the right hand side vector due to volume sources.
279  */
280  template<unsigned int dim>
281  void set_sources();
282 
283  /**
284  * @brief Assembles the fluxes on the boundary.
285  */
286  template<unsigned int dim>
288 
289  /**
290  * @brief Assembles the fluxes between elements of the same dimension.
291  */
292  template<unsigned int dim>
294 
295  /**
296  * @brief Assembles the fluxes between elements of different dimensions.
297  */
298  template<unsigned int dim>
300 
301 
302  /**
303  * @brief Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
304  *
305  * The routine just calls templated method set_boundary_condition() for each space dimension.
306  */
308 
309  /**
310  * @brief Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions
311  * for a given space dimension.
312  */
313  template<unsigned int dim>
315 
316  /**
317  * @brief Calculates the velocity field on a given @p dim dimensional cell.
318  *
319  * @param cell The cell.
320  * @param velocity The computed velocity field (at quadrature points).
321  * @param fv The FEValues class providing the quadrature points
322  * and the shape functions for velocity.
323  */
324  template<unsigned int dim>
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 up some parameters of the DG method for two sides of an edge.
344  *
345  * @param edg The edge.
346  * @param s1 Side 1.
347  * @param s2 Side 2.
348  * @param K_size Size of vector of tensors K.
349  * @param K1 Dispersivity tensors on side s1 (in quadrature points).
350  * @param K2 Dispersivity tensors on side s2 (in quadrature points).
351  * @param normal_vector Normal vector to side 0 of the neighbour
352  * (assumed constant along the side).
353  * @param alpha1, alpha2 Penalty parameter that influences the continuity
354  * of the solution (large value=more continuity).
355  * @param gamma Computed penalty parameters.
356  * @param omega Computed weights.
357  * @param transport_flux Computed flux from side s1 to side s2.
358  */
359  void set_DG_parameters_edge(const Edge &edg,
360  const int s1,
361  const int s2,
362  const int K_size,
363  const std::vector<arma::mat33> &K1,
364  const std::vector<arma::mat33> &K2,
365  const std::vector<double> &fluxes,
366  const arma::vec3 &normal_vector,
367  const double alpha1,
368  const double alpha2,
369  double &gamma,
370  double *omega,
371  double &transport_flux);
372 
373  /**
374  * @brief Sets up parameters of the DG method on a given boundary edge.
375  *
376  * Assumption is that the edge consists of only 1 side.
377  * @param side The boundary side.
378  * @param K_size Size of vector of tensors K.
379  * @param K Dispersivity tensor.
380  * @param ad_vector Advection vector.
381  * @param normal_vector Normal vector (assumed constant along the edge).
382  * @param alpha Penalty parameter that influences the continuity
383  * of the solution (large value=more continuity).
384  * @param gamma Computed penalty parameters.
385  */
386  void set_DG_parameters_boundary(const SideIter side,
387  const int K_size,
388  const std::vector<arma::mat33> &K,
389  const double flux,
390  const arma::vec3 &normal_vector,
391  const double alpha,
392  double &gamma);
393 
394 
395  /**
396  * @brief Sets the initial condition.
397  */
398  void set_initial_condition();
399 
400  /**
401  * @brief Assembles the auxiliary linear system to calculate the initial solution
402  * as L^2-projection of the prescribed initial condition.
403  */
404  template<unsigned int dim>
406 
407  /**
408  * @brief Calculates flux through boundary of each region.
409  *
410  * This actually calls calc_fluxes<dim>() for each space dimension.
411  * @param bcd_balance Total fluxes.
412  * @param bcd_plus_balance Incoming fluxes.
413  * @param bcd_minus_balance Outgoing fluxes.
414  */
415  void calc_fluxes(vector<vector<double> > &bcd_balance, vector<vector<double> > &bcd_plus_balance, vector<vector<double> > &bcd_minus_balance);
416 
417  /**
418  * @brief Calculates flux through boundary of each region of specific dimension.
419  * @param bcd_balance Total fluxes.
420  * @param bcd_plus_balance Incoming fluxes.
421  * @param bcd_minus_balance Outgoing fluxes.
422  */
423  template<unsigned int dim>
424  void calc_fluxes(vector<vector<double> > &bcd_balance, vector<vector<double> > &bcd_plus_balance, vector<vector<double> > &bcd_minus_balance);
425 
426  /**
427  * @brief Calculates volume sources for each region.
428  *
429  * This method actually calls calc_elem_sources<dim>() for each space dimension.
430  * @param mass Vector of substance mass per region.
431  * @param src_balance Vector of sources per region.
432  */
433  void calc_elem_sources(vector<vector<double> > &mass, vector< vector<double> > &src_balance);
434 
435  /**
436  * @brief Calculates volume sources for each region of specific dimension.
437  * @param mass Vector of substance mass per region.
438  * @param src_balance Vector of sources per region.
439  */
440  template<unsigned int dim>
441  void calc_elem_sources(vector<vector<double> > &mass, vector< vector<double> > &src_balance);
442 
443 
444 
445  /// @name Physical parameters
446  // @{
447 
448  /// Field data for model parameters.
450 
451  // @}
452 
453 
454  /// @name Parameters of the numerical method
455  // @{
456 
457  /// Finite element objects
459 
460  /// Penalty parameters.
462 
463  /// DG variant ((non-)symmetric/incomplete
465 
466  /// Polynomial order of finite elements.
467  unsigned int dg_order;
468 
469  // @}
470 
471 
472 
473  /// @name Solution of algebraic system
474  // @{
475 
476  /// Vector of right hand side.
477  Vec *rhs;
478 
479  /// The stiffness matrix.
481 
482  /// The mass matrix.
484 
485  /// Mass from previous time instant (necessary when coefficients of mass matrix change in time).
486  Vec *mass_vec;
487 
488  /// Linear algebra system for the transport equation.
490 
491  /// Linear algebra system for the time derivative (actually it is used only for handling the matrix structures).
493 
494  // @}
495 
496 
497  /// @name Output to file
498  // @{
499 
500  /// Array for storing the output solution data.
502 
503  /// Vector of solution data.
505 
506  /// Record with output specification.
508 
510 
511 
512  // @}
513 
514 
515  /// @name Auxiliary fields used during assembly
516  // @{
517 
518  /// Mass matrix coefficients.
520  /// Advection coefficients.
522  /// Diffusion coefficients.
524  /// Advection coefficients on edges.
526  /// Diffusion coefficients on edges.
528  /// List of indices used to call balance methods for a set of quantities.
530 
531  // @}
532 
533 
534 
535 
536  /// @name Other
537  // @{
538 
539  /// Indicates whether matrices have been preallocated.
541 
542  // @}
543 };
544 
545 
546 
547 
548 
549 
550 #endif /* TRANSPORT_DG_HH_ */
static Input::Type::Selection output_selection
void assemble_fluxes_element_element()
Assembles the fluxes between elements of the same dimension.
void set_sources()
Assembles the right hand side due to volume sources.
void set_boundary_conditions()
Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions.
FiniteElement< dim, 3 > * fe()
vector< double > mm_coef
Mass matrix coefficients.
void assemble_fluxes_element_side()
Assembles the fluxes between elements of different dimensions.
DOFHandlerMultiDim * dh_
Object for distribution of dofs.
Definition: transport_dg.hh:97
Input::Record output_rec
Record with output specification.
Transport with dispersion implemented using discontinuous Galerkin method.
Field< 3, FieldValue< 3 >::Integer > region_id
Quadrature< 1 > * q1_
Definition: transport_dg.hh:86
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.
void assemble_mass_matrix()
Assembles the mass matrix.
FiniteElement< 2, 3 > * fe_rt2_
Definition: transport_dg.hh:81
OutputTime * output_stream
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:52
void update_solution() override
Computes the solution in one time instant.
FiniteElement< 2, 3 > * fe2_
Definition: transport_dg.hh:76
int dg_variant
DG variant ((non-)symmetric/incomplete.
void prepare_initial_condition()
Assembles the auxiliary linear system to calculate the initial solution as L^2-projection of the pres...
Definition: mesh.h:109
void set_initial_condition()
Sets the initial condition.
FiniteElement< 3, 3 > * fe3_
Definition: transport_dg.hh:77
Mapping< 0, 3 > * map0_
Auxiliary mappings of reference elements.
Definition: transport_dg.hh:91
vector< vector< arma::vec3 > > ad_coef
Advection coefficients.
Definition: edges.h:38
void assemble_stiffness_matrix()
Assembles the stiffness matrix.
vector< vector< vector< arma::mat33 > > > dif_coef_edg
Diffusion coefficients on edges.
Quadrature< 2 > * q2_
Definition: transport_dg.hh:87
Specification of transport model interface.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:42
void assemble_fluxes_boundary()
Assembles the fluxes on the boundary.
TransportDG(Mesh &init_mesh, const Input::Record &in_rec)
Constructor.
vector< double * > output_solution
Array for storing the output solution data.
vector< vector< vector< arma::vec3 > > > ad_coef_edg
Advection coefficients on edges.
void output_data()
Postprocesses the solution and writes to output file.
void preallocate()
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
vector< vector< arma::mat33 > > dif_coef
Diffusion coefficients.
FEObjects * feo
Finite element objects.
FiniteElement< 1, 3 > * fe_rt1_
Finite elements for the water velocity field.
Definition: transport_dg.hh:80
Abstract class for the mapping between reference and actual cell.
Definition: fe_values.hh:44
FiniteElement< 1, 3 > * fe1_
Finite elements for the solution of the advection-diffusion equation.
Definition: transport_dg.hh:75
Accessor to the data with type Type::Record.
Definition: accessors.hh:327
Quadrature< 3 > * q3_
Definition: transport_dg.hh:88
Mat mass_matrix
The mass matrix.
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:251
Mat * stiffness_matrix
The stiffness matrix.
Field< 3, FieldValue< 3 >::Vector > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
FiniteElement< dim, 3 > * fe_rt()
unsigned int dg_order
Polynomial order of finite elements.
void output_vector_gather()
Mapping< 2, 3 > * map2_
Definition: transport_dg.hh:93
The class for outputting data during time.
Definition: output_time.hh:32
Vec * rhs
Vector of right hand side.
static Input::Type::Record input_type
Declare input record type for the equation TransportDG.
BCField< 3, FieldValue< 3 >::Vector > bc_flux
Flux in Neumann or Robin b.c.
static Input::Type::Selection bc_type_selection
DOFHandlerMultiDim * dh()
LinSys * ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
Model::ModelEqData & data()
FiniteElement< 3, 3 > * fe_rt3_
Definition: transport_dg.hh:82
EqData data_
Field data for model parameters.
bool allocation_done
Indicates whether matrices have been preallocated.
std::vector< std::vector< double > > gamma
Penalty parameters.
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.
BCField< 3, FieldValue< 3 >::EnumVector > bc_type
Type of boundary condition (see also BC_Type)
Vec * mass_vec
Mass from previous time instant (necessary when coefficients of mass matrix change in time)...
Mapping< 1, 3 > * map1_
Definition: transport_dg.hh:92
Abstract linear system class.
virtual EqData * get_data()
Getter for field data.
void calc_elem_sources(vector< vector< double > > &mass, vector< vector< double > > &src_balance)
Calculates volume sources for each region.
vector< Vec > output_vec
Vector of solution data.
Quadrature< dim > * q()
virtual void set_velocity_field(const MH_DofHandler &dh)
Updates the velocity field which determines some coefficients of the transport equation.
Record type proxy class.
Definition: type_record.hh:169
void assemble_volume_integrals()
Assembles the volume integrals into the stiffness matrix.
Mapping< 3, 3 > * map3_
Definition: transport_dg.hh:94
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: dofhandler.hh:40
Base class for FEValues and FESideValues.
Definition: fe_values.hh:43
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.
Mapping< dim, 3 > * mapping()
static Input::Type::Selection dg_variant_selection_input_type
Input type for the DG variant selection.
~TransportDG()
Destructor.
LinSys ** ls
Linear algebra system for the transport equation.
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.
Field< 3, FieldValue< 3 >::Vector > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
Template for classes storing finite set of named values.
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.
Quadrature< 0 > * q0_
Quadratures used in assembling methods.
Definition: transport_dg.hh:85
BCField< 3, FieldValue< 3 >::Vector > bc_robin_sigma
Transition coefficient in Robin b.c.
void zero_time_step() override
Initialize solution in the zero time.
FEObjects(Mesh *mesh_, unsigned int fe_order)
Definition: transport_dg.cc:96