Flow123d  JS_before_hm-1626-gde32303
transport.h
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.h
15  * @brief
16  * @todo
17  * - remove transport_sources
18  * - in create_transport_matric_mpi, there there is condition edge_flow > ZERO
19  * this makes matrix sparser, but can lead to elements without outflow and other problems
20  * when there are big differences in fluxes, more over it doesn't work if overall flow is very small
21  */
22 
23 #ifndef TRANSPORT_H_
24 #define TRANSPORT_H_
25 
26 
27 #include <memory> // for shared_ptr
28 #include <vector> // for vector
29 #include <petscmat.h>
30 #include "fem/fe_values.hh" // for FEValues
31 #include "fields/field.hh" // for Field
32 #include "fields/bc_multi_field.hh"
33 #include "fields/field_values.hh"
34 #include "fields/multi_field.hh"
35 #include "la/vector_mpi.hh"
37 #include "input/type_base.hh" // for Array
38 #include "input/type_generic.hh" // for Instance
39 #include "input/accessors.hh"
40 #include "system/index_types.hh"
41 #include "mesh/region.hh" // for RegionSet
42 #include "petscvec.h" // for Vec, _p_Vec
43 #include "tools/time_marks.hh" // for TimeMark, TimeM...
44 #include "transport/substance.hh" // for SubstanceList
47 #include "la/vector_mpi.hh"
48 
49 class OutputTime;
50 class Mesh;
51 class Distribution;
52 class Balance;
53 namespace Input {
54  namespace Type {
55  class Record;
56  class Selection;
57  }
58 }
59 
60 
61 //=============================================================================
62 // TRANSPORT
63 //=============================================================================
64 
65 /**
66  * TODO:
67  * - doxy documentation
68  * - make separate method for changing time step (rescaling only), reassembly matrix only when data are changed
69  *
70  * - needs support from EqData to determine next change of the data for 1) transport matrix 2) source vector
71  * this allows us precisely choose interval where to fix timestep
72  * : field->next_change_time() - returns time jump or actual time in case of time dep. field
73  */
74 
75 
76 
77 /**
78  * Auxiliary container class for Finite element and related objects of all dimensions.
79  * Its purpose is to provide templated access to these objects, applicable in
80  * the assembling methods.
81  */
83 public:
84 
87 
88  template<unsigned int dim>
89  inline FiniteElement<dim> *fe();
90 
91  inline Quadrature &q(unsigned int dim);
92 
93  inline FEValues<3> &fe_values(unsigned int dim)
94  {
95  ASSERT_DBG( dim >= 1 && dim <= 3 );
96  return fe_values_[dim-1];
97  }
98 
99 private:
100 
101  /// Finite elements for the solution of the advection-diffusion equation.
106 
107  /// Quadratures used in assembling methods.
109 
110  /// FESideValues objects for side flux calculating.
111  FEValues<3> fe_values_[3];
112 };
113 
114 
115 /**
116  * Class that implements explicit finite volumes scheme with upwind. The timestep is given by CFL condition.
117  *
118  */
120 public:
121  class EqData : public TransportEqFields {
122  public:
123 
124  EqData();
125  virtual ~EqData() {};
126 
127  /// Override generic method in order to allow specification of the boundary conditions through the old bcd files.
128  RegionSet read_boundary_list_item(Input::Record rec);
129 
130  /**
131  * Boundary conditions (Dirichlet) for concentrations.
132  * They are applied only on water inflow part of the boundary.
133  */
135 
136  /// Initial concentrations.
138 
141 
142  MultiField<3, FieldValue<3>::Scalar> conc_mobile; ///< Calculated concentrations in the mobile zone.
143  FieldFEScalarVec conc_mobile_fe; ///< Underlaying FieldFE for each substance of conc_mobile.
144 
145  /// Fields indended for output, i.e. all input fields plus those representing solution.
147  };
148 
149 
151 
152  static const Input::Type::Record & get_input_type();
153 
154  static const IT::Selection & get_output_selection();
155 
156  /**
157  * Constructor.
158  */
159  ConvectionTransport(Mesh &init_mesh, const Input::Record in_rec);
160  /**
161  * TODO: destructor
162  */
163  virtual ~ConvectionTransport();
164 
165  void initialize() override;
166 
167  /**
168  * Initialize solution at zero time.
169  */
170  void zero_time_step() override;
171  /**
172  * Evaluates CFL condition.
173  * Assembles the transport matrix and vector (including sources, bc terms).
174  * @param time_constraint is the value CFL constraint (return parameter)
175  * @return true if CFL is changed since previous step, false otherwise
176  */
177  bool evaluate_time_constraint(double &time_constraint) override;
178  /**
179  * Calculates one time step of explicit transport.
180  */
181  void update_solution() override;
182 
183  /** Compute P0 interpolation of the solution (used reaction term).
184  * Empty - solution is already P0 interpolation.
185  */
186  void compute_p0_interpolation() override {};
187 
188  /// Not used in this class.
189  void update_after_reactions(bool) override {};
190 
191  /**
192  * Set time interval which is considered as one time step by TransportOperatorSplitting.
193  * In particular the velocity field doesn't change over this interval.
194  *
195  * Dependencies:
196  *
197  * velocity, porosity -> matrix, source_vector
198  * matrix -> time_step
199  *
200  * data_read_times -> time_step (not necessary if we won't stick to jump times)
201  * data -> source_vector
202  * time_step -> scaling
203  *
204  *
205  *
206  */
207  void set_target_time(double target_time) override;
208 
209  /**
210  * Use Balance object from upstream equation (e.g. in various couplings) instead of own instance.
211  */
212  void set_balance_object(std::shared_ptr<Balance> balance) override;
213 
215  { return subst_idx; }
216 
217  /**
218  * @brief Write computed fields.
219  */
220  virtual void output_data() override;
221 
222  void set_output_stream(std::shared_ptr<OutputTime> stream) override
223  { output_stream_ = stream; }
224 
225 
226  /**
227  * Getters.
228  */
229 
230  /// Getter for P0 interpolation by FieldFE.
231  FieldFEScalarVec& get_p0_interpolation() override;
232 
233  Vec get_component_vec(unsigned int sbi) override;
234 
235  void get_par_info(LongIdx * &el_4_loc, Distribution * &el_ds) override;
236 
237  LongIdx *get_row_4_el() override;
238 
239  /// Returns number of transported substances.
240  inline unsigned int n_substances() override
241  { return substances_.size(); }
242 
243  /// Returns reference to the vector of substance names.
244  inline SubstanceList &substances() override
245  { return substances_; }
246 
247 private:
248 
249  /**
250  * Assembly convection term part of the matrix and boundary matrix for application of boundary conditions.
251  *
252  * Discretization of the convection term use explicit time scheme and finite volumes with full upwinding.
253  * We count on with exchange between dimensions and mixing on edges where more then two elements connect (can happen for 2D and 1D elements in
254  * 3D embedding space)
255  *
256  * In order to get multiplication matrix for explicit transport one have to scale the convection part by the acctual time step and
257  * add time term, i. e. unit matrix (see. transport_matrix_step_mpi)
258  *
259  * Updates CFL time step constrain.
260  */
261  void create_transport_matrix_mpi();
262  void create_mass_matrix();
263 
264  void make_transport_partitioning(); //
265  void set_initial_condition();
266  void read_concentration_sources();
267 
268  /** @brief Assembles concentration sources for each substance and set boundary conditions.
269  * note: the source of concentration is multiplied by time interval (gives the mass, not the flow like before)
270  */
271  void conc_sources_bdr_conditions();
272 
273  /**
274  * Finish explicit transport matrix (time step scaling)
275  */
276  void transport_matrix_step_mpi(double time_step); //
277 
278  void alloc_transport_vectors();
279  void alloc_transport_structs_mpi();
280 
281  /**
282  * @brief Wrapper of following method, call side_flux with correct template parameter.
283  */
284  double side_flux(const DHCellSide &cell_side);
285 
286  /**
287  * @brief Calculate flux on side of given element specified by dimension.
288  */
289  template<unsigned int dim>
290  double calculate_side_flux(const DHCellSide &cell);
291 
292 
293 
294  /// Registrar of class to factory
295  static const int registrar;
296 
297  /**
298  * Parameters of the equation, some are shared with other implementations since EqData is derived from TransportBase::TransportEqFields
299  */
301 
302  //@{
303  /**
304  * Flag indicates the state of object (transport matrix or source or boundary term).
305  * If false, the object is freshly assembled and not rescaled.
306  * If true, the object is scaled (not necessarily with the current time step).
307  */
308  bool is_convection_matrix_scaled, is_src_term_scaled, is_bc_term_scaled;
309 
310  /// Flag indicates that porosity or cross_section changed during last time.
312  //@}
313 
315 
316 
317  TimeMark::Type target_mark_type; ///< TimeMark type for time marks denoting end of every time interval where transport matrix remains constant.
318  double cfl_max_step; ///< Time step constraint coming from CFL condition.
319 
320  Vec vcfl_flow_, ///< Parallel vector for flow contribution to CFL condition.
321  vcfl_source_; ///< Parallel vector for source term contribution to CFL condition.
322  double *cfl_flow_, *cfl_source_;
323 
324 
325  VecScatter vconc_out_scatter;
326  Mat tm; // PETSc transport matrix
327  Vec mass_diag; // diagonal entries in pass matrix (cross_section * porosity)
328  Vec vpmass_diag; // diagonal entries in mass matrix from last time (cross_section * porosity)
329  Vec *v_tm_diag; // additions to PETSC transport matrix on the diagonal - from sources (for each substance)
330  double **tm_diag;
331 
332  /// Time when the transport matrix was created.
333  /// TODO: when we have our own classes for LA objects, we can use lazy dependence to check
334  /// necessity for matrix update
336  double transport_bc_time; ///< Time of the last update of the boundary condition terms.
337 
338  ///
339  Vec *vpconc; // previous concentration vector
340  Vec *bcvcorr; // boundary condition correction vector
342  double **cumulative_corr;
343 
344  /// Record with input specification.
346 
347  std::shared_ptr<OutputTime> output_stream_;
348 
349 
353 
354  /// Transported substances.
356 
357  /**
358  * Temporary solution how to pass velocity field form the flow model.
359  * TODO: introduce FieldDiscrete -containing true DOFHandler and data vector and pass such object together with other
360  * data. Possibly make more general set_data method, allowing setting data given by name. needs support from EqDataBase.
361  */
362  std::shared_ptr<DOFHandlerMultiDim> dh_;
363 
364  /// List of indices used to call balance methods for a set of quantities.
366 
367  /// Finite element objects
369 
371 };
372 #endif /* TRANSPORT_H_ */
bool is_mass_diag_changed
Flag indicates that porosity or cross_section changed during last time.
Definition: transport.h:311
const vector< unsigned int > & get_subst_idx() override
Return substance indices used in balance.
Definition: transport.h:214
double transport_matrix_time
Definition: transport.h:335
QGauss::array q_
Quadratures used in assembling methods.
Definition: transport.h:108
unsigned int n_substances() override
Returns number of transported substances.
Definition: transport.h:240
Abstract linear system class.
Definition: balance.hh:40
MultiField< 3, FieldValue< 3 >::Scalar > conc_mobile
Calculated concentrations in the mobile zone.
Definition: transport.h:142
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:95
std::shared_ptr< OutputTime > output_stream_
Definition: transport.h:347
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:317
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: transport.h:139
double transport_bc_time
Time of the last update of the boundary condition terms.
Definition: transport.h:336
Coupling of a transport model with a reaction model by operator splitting.
Definition: mesh.h:77
EquationOutput output_fields
Fields indended for output, i.e. all input fields plus those representing solution.
Definition: transport.h:146
LongIdx * el_4_loc
Definition: transport.h:351
Class FEValues calculates finite element data on the actual cells such as shape function values...
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_conc
Definition: transport.h:134
MultiField< 3, FieldValue< 3 >::Scalar > init_conc
Initial concentrations.
Definition: transport.h:137
const Input::Record input_rec
Record with input specification.
Definition: transport.h:345
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
SubstanceList & substances() override
Returns reference to the vector of substance names.
Definition: transport.h:244
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: transport.h:140
FieldFEScalarVec conc_mobile_fe
Underlaying FieldFE for each substance of conc_mobile.
Definition: transport.h:143
vector< VectorMPI > corr_vec
Definition: transport.h:314
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
static const int registrar
Registrar of class to factory.
Definition: transport.h:295
SubstanceList substances_
Transported substances.
Definition: transport.h:355
FEValues< 3 > & fe_values(unsigned int dim)
Definition: transport.h:93
std::shared_ptr< DOFHandlerMultiDim > dh_
Definition: transport.h:362
The class for outputting data during time.
Definition: output_time.hh:51
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
Definition: transport.h:365
void update_after_reactions(bool) override
Not used in this class.
Definition: transport.h:189
void compute_p0_interpolation() override
Definition: transport.h:186
std::array< QGauss, 4 > array
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
double * cfl_source_
Definition: transport.h:322
void set_output_stream(std::shared_ptr< OutputTime > stream) override
Setter for output stream.
Definition: transport.h:222
VecScatter vconc_out_scatter
Definition: transport.h:325
FiniteElement< 1 > * fe1_
Definition: transport.h:103
double ** tm_diag
Definition: transport.h:330
ConcentrationTransportBase FactoryBaseType
Definition: transport.h:150
Classes for storing substance data.
#define ASSERT_DBG(expr)
Vec vcfl_source_
Parallel vector for source term contribution to CFL condition.
Definition: transport.h:320
Definitions of particular quadrature rules on simplices.
double ** cumulative_corr
Definition: transport.h:342
Record type proxy class.
Definition: type_record.hh:182
FETransportObjects feo_
Finite element objects.
Definition: transport.h:368
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...
FiniteElement< 3 > * fe3_
Definition: transport.h:105
FiniteElement< 0 > * fe0_
Finite elements for the solution of the advection-diffusion equation.
Definition: transport.h:102
Distribution * el_ds
Definition: transport.h:352
Template for classes storing finite set of named values.
Side accessor allows to iterate over sides of DOF handler cell.
LongIdx * row_4_el
Definition: transport.h:350
double cfl_max_step
Time step constraint coming from CFL condition.
Definition: transport.h:318
FiniteElement< 2 > * fe2_
Definition: transport.h:104