Flow123d  JS_before_hm-1-ge159291
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 #include <boost/exception/info.hpp> // for operator<<, err...
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 "mesh/long_idx.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
45 
46 class OutputTime;
47 class Mesh;
48 class Distribution;
49 class Balance;
50 class MH_DofHandler;
51 namespace Input {
52  namespace Type {
53  class Record;
54  class Selection;
55  }
56 }
57 template <int spacedim, class Value> class FieldFE;
58 
59 
60 //=============================================================================
61 // TRANSPORT
62 //=============================================================================
63 
64 /**
65  * TODO:
66  * - doxy documentation
67  * - make separate method for changing time step (rescaling only), reassembly matrix only when data are changed
68  *
69  * - needs support from EqData to determine next change of the data for 1) transport matrix 2) source vector
70  * this allows us precisely choose interval where to fix timestep
71  * : field->next_change_time() - returns time jump or actual time in case of time dep. field
72  */
73 
74 
75 
76 /**
77  * Class that implements explicit finite volumes scheme with upwind. The timestep is given by CFL condition.
78  *
79  */
81 public:
82 
83  class EqData : public TransportEqData {
84  public:
85 
86  EqData();
87  virtual ~EqData() {};
88 
89  /// Override generic method in order to allow specification of the boundary conditions through the old bcd files.
90  RegionSet read_boundary_list_item(Input::Record rec);
91 
92  /**
93  * Boundary conditions (Dirichlet) for concentrations.
94  * They are applied only on water inflow part of the boundary.
95  */
97 
98  /// Initial concentrations.
100 
103  MultiField<3, FieldValue<3>::Scalar> conc_mobile; ///< Calculated concentrations in the mobile zone.
104 
105 
106  /// Fields indended for output, i.e. all input fields plus those representing solution.
108  };
109 
110 
112 
113  static const Input::Type::Record & get_input_type();
114 
115  static const IT::Selection & get_output_selection();
116 
117  /**
118  * Constructor.
119  */
120  ConvectionTransport(Mesh &init_mesh, const Input::Record in_rec);
121  /**
122  * TODO: destructor
123  */
124  virtual ~ConvectionTransport();
125 
126  void initialize() override;
127 
128  /**
129  * Initialize solution at zero time.
130  */
131  void zero_time_step() override;
132  /**
133  * Evaluates CFL condition.
134  * Assembles the transport matrix and vector (including sources, bc terms).
135  * @param time_constraint is the value CFL constraint (return parameter)
136  * @return true if CFL is changed since previous step, false otherwise
137  */
138  bool evaluate_time_constraint(double &time_constraint) override;
139  /**
140  * Calculates one time step of explicit transport.
141  */
142  void update_solution() override;
143 
144  void calculate_concentration_matrix() override {};
145 
146  /// Not used in this class.
147  void update_after_reactions(bool) override {};
148 
149  /**
150  * Set time interval which is considered as one time step by TransportOperatorSplitting.
151  * In particular the velocity field doesn't change over this interval.
152  *
153  * Dependencies:
154  *
155  * velocity, porosity -> matrix, source_vector
156  * matrix -> time_step
157  *
158  * data_read_times -> time_step (not necessary if we won't stick to jump times)
159  * data -> source_vector
160  * time_step -> scaling
161  *
162  *
163  *
164  */
165  void set_target_time(double target_time) override;
166 
167  /**
168  * Use Balance object from upstream equation (e.g. in various couplings) instead of own instance.
169  */
170  void set_balance_object(std::shared_ptr<Balance> balance) override;
171 
173  { return subst_idx; }
174 
175  /**
176  * @brief Write computed fields.
177  */
178  virtual void output_data() override;
179 
180  inline void set_velocity_field(const MH_DofHandler &dh) override
181  { mh_dh=&dh; }
182 
183  void set_output_stream(std::shared_ptr<OutputTime> stream) override
184  { output_stream_ = stream; }
185 
186 
187  /**
188  * Getters.
189  */
190  inline std::shared_ptr<OutputTime> output_stream() override
191  { return output_stream_; }
192 
193  double **get_concentration_matrix() override;
194 
195  const Vec &get_solution(unsigned int sbi) override
196  { return vconc[sbi]; }
197 
198  void get_par_info(LongIdx * &el_4_loc, Distribution * &el_ds) override;
199 
200  LongIdx *get_row_4_el() override;
201 
202  /// Returns number of transported substances.
203  inline unsigned int n_substances() override
204  { return substances_.size(); }
205 
206  /// Returns reference to the vector of substance names.
207  inline SubstanceList &substances() override
208  { return substances_; }
209 
210 private:
211 
212  /**
213  * Assembly convection term part of the matrix and boundary matrix for application of boundary conditions.
214  *
215  * Discretization of the convection term use explicit time scheme and finite volumes with full upwinding.
216  * 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
217  * 3D embedding space)
218  *
219  * In order to get multiplication matrix for explicit transport one have to scale the convection part by the acctual time step and
220  * add time term, i. e. unit matrix (see. transport_matrix_step_mpi)
221  *
222  * Updates CFL time step constrain.
223  */
224  void create_transport_matrix_mpi();
225  void create_mass_matrix();
226 
227  void make_transport_partitioning(); //
228  void set_initial_condition();
229  void read_concentration_sources();
230  void set_boundary_conditions();
231 
232  /** @brief Assembles concentration sources for each substance.
233  * note: the source of concentration is multiplied by time interval (gives the mass, not the flow like before)
234  */
235  void compute_concentration_sources();
236 
237  /**
238  * Finish explicit transport matrix (time step scaling)
239  */
240  void transport_matrix_step_mpi(double time_step); //
241 
242  void alloc_transport_vectors();
243  void alloc_transport_structs_mpi();
244 
245  /**
246  * Communicate parallel concentration vectors into sequential output vector.
247  */
248  void output_vector_gather();
249 
250 
251 
252  /// Registrar of class to factory
253  static const int registrar;
254 
255  /**
256  * Parameters of the equation, some are shared with other implementations since EqData is derived from TransportBase::TransportEqData
257  */
259 
260  //@{
261  /**
262  * Flag indicates the state of object (transport matrix or source or boundary term).
263  * If false, the object is freshly assembled and not rescaled.
264  * If true, the object is scaled (not necessarily with the current time step).
265  */
266  bool is_convection_matrix_scaled, is_src_term_scaled, is_bc_term_scaled;
267 
268  /// Flag indicates that porosity or cross_section changed during last time.
270  //@}
271 
272  double **sources_corr;
274 
275 
276  TimeMark::Type target_mark_type; ///< TimeMark type for time marks denoting end of every time interval where transport matrix remains constant.
277  double cfl_max_step; ///< Time step constraint coming from CFL condition.
278 
279  Vec vcfl_flow_, ///< Parallel vector for flow contribution to CFL condition.
280  vcfl_source_; ///< Parallel vector for source term contribution to CFL condition.
281  double *cfl_flow_, *cfl_source_;
282 
283 
284  VecScatter vconc_out_scatter;
285  Mat tm; // PETSc transport matrix
286  Vec mass_diag; // diagonal entries in pass matrix (cross_section * porosity)
287  Vec vpmass_diag; // diagonal entries in mass matrix from last time (cross_section * porosity)
288  Vec *v_tm_diag; // additions to PETSC transport matrix on the diagonal - from sources (for each substance)
289  double **tm_diag;
290 
291  /// Time when the transport matrix was created.
292  /// TODO: when we have our own classes for LA objects, we can use lazy dependence to check
293  /// necessity for matrix update
295  double transport_bc_time; ///< Time of the last update of the boundary condition terms.
296 
297  /// Concentration vectors for mobile phase.
298  Vec *vconc; // concentration vector
299  /// Concentrations for phase, substance, element
300  double **conc;
301 
302  ///
303  Vec *vpconc; // previous concentration vector
304  Vec *bcvcorr; // boundary condition correction vector
306  double **cumulative_corr;
307 
309 
310  // Temporary objects holding pointers to appropriate FieldFE
311  // TODO remove after final fix of equations
312  /// Fields correspond with \p out_conc.
314 
315  /// Record with input specification.
317 
318  std::shared_ptr<OutputTime> output_stream_;
319 
320 
324 
325  /// Transported substances.
327 
328  /**
329  * Temporary solution how to pass velocity field form the flow model.
330  * TODO: introduce FieldDiscrete -containing true DOFHandler and data vector and pass such object together with other
331  * data. Possibly make more general set_data method, allowing setting data given by name. needs support from EqDataBase.
332  */
334  std::shared_ptr<DOFHandlerMultiDim> dh_;
335 
336  /// List of indices used to call balance methods for a set of quantities.
338 
340 };
341 #endif /* TRANSPORT_H_ */
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
bool is_mass_diag_changed
Flag indicates that porosity or cross_section changed during last time.
Definition: transport.h:269
const vector< unsigned int > & get_subst_idx() override
Return substance indices used in balance.
Definition: transport.h:172
double transport_matrix_time
Definition: transport.h:294
const Vec & get_solution(unsigned int sbi) override
Return PETSc vector with solution for sbi-th substance.
Definition: transport.h:195
unsigned int n_substances() override
Returns number of transported substances.
Definition: transport.h:203
Abstract linear system class.
Definition: balance.hh:37
MultiField< 3, FieldValue< 3 >::Scalar > conc_mobile
Calculated concentrations in the mobile zone.
Definition: transport.h:103
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:83
void set_velocity_field(const MH_DofHandler &dh) override
Pass velocity from flow to transport.
Definition: transport.h:180
std::shared_ptr< OutputTime > output_stream_
Definition: transport.h:318
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:276
Field< 3, FieldValue< 3 >::Scalar > region_id
Definition: transport.h:101
double transport_bc_time
Time of the last update of the boundary condition terms.
Definition: transport.h:295
Coupling of a transport model with a reaction model by operator splitting.
Definition: mesh.h:76
EquationOutput output_fields
Fields indended for output, i.e. all input fields plus those representing solution.
Definition: transport.h:107
double ** sources_corr
Definition: transport.h:272
LongIdx * el_4_loc
Definition: transport.h:322
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_conc
Definition: transport.h:96
std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > output_field_ptr
Fields correspond with out_conc.
Definition: transport.h:313
double ** conc
Concentrations for phase, substance, element.
Definition: transport.h:300
MultiField< 3, FieldValue< 3 >::Scalar > init_conc
Initial concentrations.
Definition: transport.h:99
const Input::Record input_rec
Record with input specification.
Definition: transport.h:316
SubstanceList & substances() override
Returns reference to the vector of substance names.
Definition: transport.h:207
Field< 3, FieldValue< 3 >::Scalar > subdomain
Definition: transport.h:102
const MH_DofHandler * mh_dh
Definition: transport.h:333
std::shared_ptr< OutputTime > output_stream() override
Definition: transport.h:190
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
void calculate_concentration_matrix() override
Calculate the array of concentrations per element (for reactions).
Definition: transport.h:144
static const int registrar
Registrar of class to factory.
Definition: transport.h:253
std::vector< VectorMPI > out_conc
Definition: transport.h:308
SubstanceList substances_
Transported substances.
Definition: transport.h:326
std::shared_ptr< DOFHandlerMultiDim > dh_
Definition: transport.h:334
The class for outputting data during time.
Definition: output_time.hh:50
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
Definition: transport.h:337
Vec * vconc
Concentration vectors for mobile phase.
Definition: transport.h:298
void update_after_reactions(bool) override
Not used in this class.
Definition: transport.h:147
double * cfl_source_
Definition: transport.h:281
void set_output_stream(std::shared_ptr< OutputTime > stream) override
Setter for output stream.
Definition: transport.h:183
VecScatter vconc_out_scatter
Definition: transport.h:284
double ** tm_diag
Definition: transport.h:289
ConcentrationTransportBase FactoryBaseType
Definition: transport.h:111
Classes for storing substance data.
Vec vcfl_source_
Parallel vector for source term contribution to CFL condition.
Definition: transport.h:279
double ** cumulative_corr
Definition: transport.h:306
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
Distribution * el_ds
Definition: transport.h:323
Template for classes storing finite set of named values.
LongIdx * row_4_el
Definition: transport.h:321
double cfl_max_step
Time step constraint coming from CFL condition.
Definition: transport.h:277
Definition: field.hh:56