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