Flow123d  DF_patch_fe_values-dbc06cd
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 "fields/field.hh" // for Field
31 #include "fields/bc_multi_field.hh"
32 #include "fields/field_values.hh"
33 #include "fields/multi_field.hh"
34 #include "la/vector_mpi.hh"
36 #include "input/type_base.hh" // for Array
37 #include "input/type_generic.hh" // for Instance
38 #include "input/accessors.hh"
39 #include "system/index_types.hh"
40 #include "mesh/region.hh" // for RegionSet
41 #include "petscvec.h" // for Vec, _p_Vec
42 #include "tools/time_marks.hh" // for TimeMark, TimeM...
43 #include "transport/substance.hh" // for SubstanceList
46 #include "la/vector_mpi.hh"
47 
48 class OutputTime;
49 class Mesh;
50 class Distribution;
51 class Balance;
52 namespace Input {
53  namespace Type {
54  class Record;
55  class Selection;
56  }
57 }
58 template<unsigned int dim, class EqData> class MassAssemblyConvection;
59 template<unsigned int dim, class EqData> class InitCondAssemblyConvection;
60 template<unsigned int dim, class EqData> class ConcSourcesBdrAssemblyConvection;
61 template<unsigned int dim, class EqData> class MatrixMpiAssemblyConvection;
62 template< template<IntDim...> class DimAssembly> class GenericAssembly;
63 
64 
65 //=============================================================================
66 // TRANSPORT
67 //=============================================================================
68 
69 /**
70  * TODO:
71  * - doxy documentation
72  * - make separate method for changing time step (rescaling only), reassembly matrix only when data are changed
73  *
74  * - needs support from EqData to determine next change of the data for 1) transport matrix 2) source vector
75  * this allows us precisely choose interval where to fix timestep
76  * : field->next_change_time() - returns time jump or actual time in case of time dep. field
77  */
78 
79 
80 
81 /**
82  * Class that implements explicit finite volumes scheme with upwind. The timestep is given by CFL condition.
83  *
84  */
86 public:
87  class EqFields : public TransportEqFields {
88  public:
89 
90  EqFields();
91  virtual ~EqFields() {};
92 
93  /// Calculate flux on given side point.
94  inline double side_flux(SidePoint &side_p, arma::vec3 normal_vec, double JxW) {
95  return arma::dot(this->flow_flux(side_p), normal_vec) * JxW;
96  }
97 
98  /**
99  * Boundary conditions (Dirichlet) for concentrations.
100  * They are applied only on water inflow part of the boundary.
101  */
103 
104  /// Initial concentrations.
106 
109 
110  MultiField<3, FieldValue<3>::Scalar> conc_mobile; ///< Calculated concentrations in the mobile zone.
111  FieldFEScalarVec conc_mobile_fe; ///< Underlaying FieldFE for each substance of conc_mobile.
112 
113  /// Fields indended for output, i.e. all input fields plus those representing solution.
115  };
116 
117 
118  class EqData {
119  public:
121 
122  EqData(std::shared_ptr<EqFields> eq_fields)
123  : eq_fields_(eq_fields), is_mass_diag_changed(false), cfl_source_(PETSC_COMM_WORLD), cfl_flow_(PETSC_COMM_WORLD) {}
124  virtual ~EqData() {};
125 
126  /// Returns number of transported substances.
127  inline unsigned int n_substances()
128  { return substances_.size(); }
129 
131  ASSERT_PTR(time);
132  this->time_ = time;
133  }
134 
135  void alloc_transport_structs_mpi(unsigned int lsize) {
136  this->cfl_flow_.resize(lsize);
137  this->cfl_source_.resize(lsize);
138  }
139 
140  /// Shared pointer of EqFields
141  std::shared_ptr<EqFields> eq_fields_;
142 
143  /**
144  * Temporary solution how to pass velocity field form the flow model.
145  * TODO: introduce FieldDiscrete -containing true DOFHandler and data vector and pass such object together with other
146  * data. Possibly make more general set_data method, allowing setting data given by name. needs support from EqDataBase.
147  */
148  std::shared_ptr<DOFHandlerMultiDim> dh_;
149 
150  /// object for calculation and writing the mass balance to file, shared with EquationBase.
151  std::shared_ptr<Balance> balance_;
152 
153  /// Flag indicates that porosity or cross_section changed during last time.
155 
156  /// Transported substances.
158 
159  /// List of indices used to call balance methods for a set of quantities.
161 
162  Vec mass_diag; // diagonal entries in pass matrix (cross_section * porosity)
163 
164  /// Flag indicates that sources part of equation was changed during last time.
166 
168  VectorMPI cfl_source_; ///< Parallel vector for source term contribution to CFL condition.
169  VectorMPI cfl_flow_; ///< Parallel vector for flow contribution to CFL condition.
170  Mat tm; ///< PETSc transport matrix
171  vector<VectorMPI> tm_diag; ///< additions to PETSC transport matrix on the diagonal - from sources (for each substance)
172  Vec *bcvcorr; ///< Boundary condition correction vector
173  double transport_bc_time; ///< Time of the last update of the boundary condition terms.
175 
176  /// Time when the transport matrix was created.
177  /// TODO: when we have our own classes for LA objects, we can use lazy dependence to check
178  /// necessity for matrix update
180 
181  bool is_convection_matrix_scaled; ///< Flag indicates the state of object
182 
183  /// Maximal number of edge sides (evaluate from dim 1,2,3)
184  unsigned int max_edg_sides;
185 
186  };
187 
192 
193 
195 
196  static const Input::Type::Record & get_input_type();
197 
199 
200  /**
201  * Constructor.
202  */
203  ConvectionTransport(Mesh &init_mesh, const Input::Record in_rec);
204  /**
205  * TODO: destructor
206  */
207  virtual ~ConvectionTransport();
208 
209  void initialize() override;
210 
211  /**
212  * Initialize solution at zero time.
213  */
214  void zero_time_step() override;
215  /**
216  * Evaluates CFL condition.
217  * Assembles the transport matrix and vector (including sources, bc terms).
218  * @param time_constraint is the value CFL constraint (return parameter)
219  * @return true if CFL is changed since previous step, false otherwise
220  */
221  bool evaluate_time_constraint(double &time_constraint) override;
222  /**
223  * Calculates one time step of explicit transport.
224  */
225  void update_solution() override;
226 
227  /** Compute P0 interpolation of the solution (used reaction term).
228  * Empty - solution is already P0 interpolation.
229  */
230  void compute_p0_interpolation() override {};
231 
232  /// Not used in this class.
233  void update_after_reactions(bool) override {};
234 
235  /**
236  * Set time interval which is considered as one time step by TransportOperatorSplitting.
237  * In particular the velocity field doesn't change over this interval.
238  *
239  * Dependencies:
240  *
241  * velocity, porosity -> matrix, source_vector
242  * matrix -> time_step
243  *
244  * data_read_times -> time_step (not necessary if we won't stick to jump times)
245  * data -> source_vector
246  * time_step -> scaling
247  *
248  *
249  *
250  */
251  void set_target_time(double target_time) override;
252 
253  /**
254  * Use Balance object from upstream equation (e.g. in various couplings) instead of own instance.
255  */
256  void set_balance_object(std::shared_ptr<Balance> balance) override;
257 
259  { return eq_data_->subst_idx; }
260 
261  /**
262  * @brief Write computed fields.
263  */
264  virtual void output_data() override;
265 
266  void set_output_stream(std::shared_ptr<OutputTime> stream) override
267  { output_stream_ = stream; }
268 
269 
270  /**
271  * Getters.
272  */
273 
274  /// Getter for P0 interpolation by FieldFE.
276 
277  Vec get_component_vec(unsigned int sbi) override;
278 
279  /// Returns number of transported substances.
280  inline unsigned int n_substances() override
281  { return eq_data_->substances_.size(); }
282 
283  /// Returns reference to the vector of substance names.
284  inline SubstanceList &substances() override
285  { return eq_data_->substances_; }
286 
287 private:
288 
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 EqFields is derived from TransportBase::TransportEqFields
299  */
300  std::shared_ptr<EqFields> eq_fields_;
301  std::shared_ptr<EqData> eq_data_;
302 
303  //@{
304  /**
305  * Flag indicates the state of object (transport matrix or source or boundary term).
306  * If false, the object is freshly assembled and not rescaled.
307  * If true, the object is scaled (not necessarily with the current time step).
308  */
310 
311  //@}
312 
313  TimeMark::Type target_mark_type; ///< TimeMark type for time marks denoting end of every time interval where transport matrix remains constant.
314  double cfl_max_step; ///< Time step constraint coming from CFL condition.
315 
316 
317  VecScatter vconc_out_scatter;
318  Vec vpmass_diag; // diagonal entries in mass matrix from last time (cross_section * porosity)
319 
320  ///
321  Vec *vpconc; // previous concentration vector
323 
324  /// Record with input specification.
326 
327  std::shared_ptr<OutputTime> output_stream_;
328 
329 
330  /// general assembly objects, hold assembly objects of appropriate dimension
335 
337 };
338 #endif /* TRANSPORT_H_ */
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR) only for debug mode.
Definition: asserts.hh:341
std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > FieldFEScalarVec
bool sources_changed_
Flag indicates that sources part of equation was changed during last time.
Definition: transport.h:165
bool is_convection_matrix_scaled
Flag indicates the state of object.
Definition: transport.h:181
void alloc_transport_structs_mpi(unsigned int lsize)
Definition: transport.h:135
Vec * bcvcorr
Boundary condition correction vector.
Definition: transport.h:172
VectorMPI cfl_source_
Parallel vector for source term contribution to CFL condition.
Definition: transport.h:168
std::shared_ptr< Balance > balance_
object for calculation and writing the mass balance to file, shared with EquationBase.
Definition: transport.h:151
std::shared_ptr< EqFields > eq_fields_
Shared pointer of EqFields.
Definition: transport.h:141
ConvectionTransport::EqFields EqFields
Definition: transport.h:120
VectorMPI cfl_flow_
Parallel vector for flow contribution to CFL condition.
Definition: transport.h:169
unsigned int max_edg_sides
Maximal number of edge sides (evaluate from dim 1,2,3)
Definition: transport.h:184
double transport_bc_time
Time of the last update of the boundary condition terms.
Definition: transport.h:173
EqData(std::shared_ptr< EqFields > eq_fields)
Definition: transport.h:122
void set_time_governor(TimeGovernor *time)
Definition: transport.h:130
bool is_mass_diag_changed
Flag indicates that porosity or cross_section changed during last time.
Definition: transport.h:154
unsigned int n_substances()
Returns number of transported substances.
Definition: transport.h:127
vector< VectorMPI > tm_diag
additions to PETSC transport matrix on the diagonal - from sources (for each substance)
Definition: transport.h:171
Mat tm
PETSc transport matrix.
Definition: transport.h:170
vector< VectorMPI > corr_vec
Definition: transport.h:167
SubstanceList substances_
Transported substances.
Definition: transport.h:157
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
Definition: transport.h:160
std::shared_ptr< DOFHandlerMultiDim > dh_
Definition: transport.h:148
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: transport.h:108
double side_flux(SidePoint &side_p, arma::vec3 normal_vec, double JxW)
Calculate flux on given side point.
Definition: transport.h:94
MultiField< 3, FieldValue< 3 >::Scalar > init_conc
Initial concentrations.
Definition: transport.h:105
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: transport.h:107
FieldFEScalarVec conc_mobile_fe
Underlaying FieldFE for each substance of conc_mobile.
Definition: transport.h:111
EquationOutput output_fields
Fields indended for output, i.e. all input fields plus those representing solution.
Definition: transport.h:114
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_conc
Definition: transport.h:102
MultiField< 3, FieldValue< 3 >::Scalar > conc_mobile
Calculated concentrations in the mobile zone.
Definition: transport.h:110
ConcentrationTransportBase FactoryBaseType
Definition: transport.h:194
void set_output_stream(std::shared_ptr< OutputTime > stream) override
Setter for output stream.
Definition: transport.h:266
SubstanceList & substances() override
Returns reference to the vector of substance names.
Definition: transport.h:284
static const IT::Selection & get_output_selection()
FieldFEScalarVec & get_p0_interpolation() override
Getter for P0 interpolation by FieldFE.
Definition: transport.cc:521
GenericAssembly< MassAssemblyConvectionDim > * mass_assembly_
general assembly objects, hold assembly objects of appropriate dimension
Definition: transport.h:331
unsigned int n_substances() override
Returns number of transported substances.
Definition: transport.h:280
VecScatter vconc_out_scatter
Definition: transport.h:317
void initialize() override
Definition: transport.cc:147
void set_target_time(double target_time) override
Definition: transport.cc:498
void update_solution() override
Definition: transport.cc:394
GenericAssembly< InitCondAssemblyConvectionDim > * init_cond_assembly_
Definition: transport.h:332
void alloc_transport_structs_mpi()
Definition: transport.cc:250
std::shared_ptr< EqFields > eq_fields_
Definition: transport.h:300
void zero_time_step() override
Definition: transport.cc:293
void alloc_transport_vectors()
Definition: transport.cc:242
void update_after_reactions(bool) override
Not used in this class.
Definition: transport.h:233
ConvectionTransport(Mesh &init_mesh, const Input::Record in_rec)
Definition: transport.cc:123
const Input::Record input_rec
Record with input specification.
Definition: transport.h:325
virtual void output_data() override
Write computed fields.
Definition: transport.cc:526
Vec get_component_vec(unsigned int sbi) override
Return PETSc vector with solution for sbi-th substance.
Definition: transport.cc:198
std::shared_ptr< OutputTime > output_stream_
Definition: transport.h:327
void set_balance_object(std::shared_ptr< Balance > balance) override
Definition: transport.cc:538
static const Input::Type::Record & get_input_type()
Definition: transport.cc:72
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:313
void compute_p0_interpolation() override
Definition: transport.h:230
GenericAssembly< ConcSourcesBdrAssemblyConvectionDim > * conc_sources_bdr_assembly_
Definition: transport.h:333
GenericAssembly< MatrixMpiAssemblyConvectionDim > * matrix_mpi_assembly_
Definition: transport.h:334
bool evaluate_time_constraint(double &time_constraint) override
Definition: transport.cc:325
double cfl_max_step
Time step constraint coming from CFL condition.
Definition: transport.h:314
static const int registrar
Registrar of class to factory.
Definition: transport.h:295
std::shared_ptr< EqData > eq_data_
Definition: transport.h:301
const vector< unsigned int > & get_subst_idx() override
Return substance indices used in balance.
Definition: transport.h:258
virtual ~ConvectionTransport()
Definition: transport.cc:205
std::shared_ptr< Balance > balance() const
Definition: equation.hh:189
TimeGovernor & time()
Definition: equation.hh:151
Class template representing a field with values dependent on: point, element, and region.
Definition: field.hh:92
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
General point a+ side_begin_ + ccessor allow iterate over quadrature points of given side defined in ...
unsigned int size() const
Definition: substance.hh:87
Basic time management functionality for unsteady (and steady) solvers (class Equation).
Field< 3, FieldValue< 3 >::VectorFixed > flow_flux
Flow flux, can be result of water flow model.
Coupling of a transport model with a reaction model by operator splitting.
void resize(unsigned int local_size)
Definition: vector_mpi.cc:50
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19
Abstract linear system class.
Definition: balance.hh:40
Definitions of particular quadrature rules on simplices.
Classes for storing substance data.