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