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