Flow123d  jenkins-Flow123d-windows32-release-multijob-51
transport.h
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief ???
27  *
28  *
29  * TODO:
30  * - remove transport_sources
31  * - in create_transport_matric_mpi, there there is condition edge_flow > ZERO
32  * this makes matrix sparser, but can lead to elements without outflow and other problems
33  * when there are big differences in fluxes, more over it doesn't work if overall flow is very small
34  *
35  */
36 
37 #ifndef TRANSPORT_H_
38 #define TRANSPORT_H_
39 
40 #include <petscmat.h>
41 #include "coupling/equation.hh"
42 #include "input/accessors.hh"
43 #include "flow/mh_dofhandler.hh"
45 
47 #include "fields/field_values.hh"
48 
49 
50 class SorptionImmob;
51 class OutputTime;
52 class Mesh;
53 class Distribution;
55 
56 //=============================================================================
57 // TRANSPORT
58 //=============================================================================
59 #define MOBILE 0
60 #define IMMOBILE 1
61 #define MOBILE_SORB 2
62 #define IMMOBILE_SORB 3
63 
64 #define MAX_PHASES 4
65 
66 /**
67  * TODO:
68  * - doxy documentation
69  * - make separate method for changing time step (rescaling only), reassembly matrix only when data are changed
70  *
71  * - needs support from EqData to determine next change of the data for 1) transport matrix 2) source vector
72  * this allows us precisely choose interval where to fix timestep
73  * : field->next_change_time() - returns time jump or actual time in case of time dep. field
74  */
75 
76 
77 
78 /**
79  * Class that implements explicit finite volumes scheme with upwind. The timestep is given by CFL condition.
80  *
81  */
83 public:
84 
86  public:
88 
90 
91  EqData();
92  virtual ~EqData() {};
93 
94  /// Override generic method in order to allow specification of the boundary conditions through the old bcd files.
96 
97  /**
98  * Boundary conditions (Dirichlet) for concentrations.
99  * They are applied only on water inflow part of the boundary.
100  */
102 
103  /// Initial concentrations.
105 
106  MultiField<3, FieldValue<3>::Scalar> conc_mobile; ///< Calculated concentrations in the mobile zone.
107 
108  /// Fields indended for output, i.e. all input fields plus those representing solution.
110  };
111 
112  /**
113  * Constructor.
114  */
115  ConvectionTransport(Mesh &init_mesh, const Input::Record &in_rec);
116  /**
117  * TODO: destructor
118  */
119  virtual ~ConvectionTransport();
120 
121  /**
122  * Initialize solution at zero time.
123  */
124  void zero_time_step() override;
125  /**
126  * Calculates one time step of explicit transport.
127  */
128  void update_solution() override;
129 
130  /**
131  * Set time interval which is considered as one time step by TransportOperatorSplitting.
132  * In particular the velocity field dosn't change over this interval.
133  *
134  * Dependencies:
135  *
136  * velocity, porosity -> matrix, source_vector
137  * matrix -> time_step
138  *
139  * data_read_times -> time_step (not necessary if we won't stick to jump times)
140  * data -> source_vector
141  * time_step -> scaling
142  *
143  *
144  *
145  */
146  void set_target_time(double target_time);
147 
148  /**
149  * Communicate parallel concentration vectors into sequential output vector.
150  */
151  void output_vector_gather();
152 
153  /**
154  * @brief Write computed fields.
155  */
156  virtual void output_data() override;
157 
158 
159  /**
160  * Getters.
161  */
162  inline EqData *get_data() { return &data_; }
163 
165 
166  double ***get_concentration_matrix();
167  void get_par_info(int * &el_4_loc, Distribution * &el_ds);
168  int *get_el_4_loc();
169  int *get_row_4_el();
170 
172 
173 private:
174 
175  /**
176  * Assembly convection term part of the matrix and boundary matrix for application of boundary conditions.
177  *
178  * Discretization of the convection term use explicit time scheme and finite volumes with full upwinding.
179  * 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
180  * 3D embedding space)
181  *
182  * In order to get multiplication matrix for explicit transport one have to scale the convection part by the acctual time step and
183  * add time term, i. e. unit matrix (see. transport_matrix_step_mpi)
184  *
185  * Updates CFL time step constrain.
186  */
188 
189  void make_transport_partitioning(); //
190  void set_initial_condition();
193 
194  //note: the source of concentration is multiplied by time interval (gives the mass, not the flow like before)
195  void compute_concentration_sources(unsigned int sbi);
196  void compute_concentration_sources_for_mass_balance(unsigned int sbi);
197 
198  /**
199  * Finish explicit transport matrix (time step scaling)
200  */
201  void transport_matrix_step_mpi(double time_step); //
202 
205 
206  /**
207  * Implements the virtual method EquationForMassBalance::calc_fluxes().
208  */
209  void calc_fluxes(vector<vector<double> > &bcd_balance, vector<vector<double> > &bcd_plus_balance, vector<vector<double> > &bcd_minus_balance);
210  /**
211  * Implements the virtual method EquationForMassBalance::calc_elem_sources().
212  */
213  void calc_elem_sources(vector<vector<double> > &mass, vector<vector<double> > &src_balance);
214 
215 
216  /**
217  * Parameters of the equation, some are shared with other implementations since EqData is derived from TransportBase::TransportEqData
218  */
220 
221  /**
222  * Indicates if we finished the matrix and add vector by scaling with timestep factor.
223  */
225 
226  //TODO: remove this and make concentration_matrix only two-dimensional
227  int sub_problem; // 0-only transport,1-transport+dual porosity,
228  // 2-transport+sorption
229  // 3-transport+dual porosity+sorption
230 
231 
232  double *sources_corr;
234 
235  ///temporary arrays to store constant values of fields over time interval
236  //(avoiding calling "field.value()" too often)
237  double **sources_density,
238  **sources_conc,
239  **sources_sigma;
240 
241  TimeMark::Type target_mark_type; ///< TimeMark type for time marks denoting end of every time interval where transport matrix remains constant.
242  double cfl_max_step;
243  // only local part
244 
245 
246 
247  VecScatter vconc_out_scatter;
248  Mat tm; // PETSc transport matrix
249 
250  /// Time when the transport matrix was created.
251  /// TODO: when we have our own classes for LA objects, we can use lazy dependence to check
252  /// necessity for matrix update
254 
255  /// Concentration vectors for mobile phase.
256  Vec *vconc; // concentration vector
257  /// Concentrations for phase, substance, element
258  double ***conc;
259 
260  ///
261  Vec *vpconc; // previous concentration vector
262  Vec *bcvcorr; // boundary condition correction vector
264  double **cumulative_corr;
265 
266  Vec *vconc_out; // concentration vector output (gathered)
267  double ***out_conc;
268 
269  /// Record with output specification.
271 
273 
274 
275  int *row_4_el;
276  int *el_4_loc;
278 
280 };
281 #endif /* TRANSPORT_H_ */
void compute_concentration_sources(unsigned int sbi)
Definition: transport.cc:370
Sorption model in immobile zone in case dual porosity is considered.
Definition: sorption.hh:97
double *** get_concentration_matrix()
Definition: transport.cc:720
Abstract base class for equation clasess.
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:51
unsigned long int Type
Definition: time_marks.hh:59
TimeIntegrationScheme time_scheme() override
Returns the time integration scheme of the equation.
Definition: transport.h:171
double * sources_corr
Definition: transport.h:232
void update_solution() override
Definition: transport.cc:484
double transport_matrix_time
Definition: transport.h:253
double *** conc
Concentrations for phase, substance, element.
Definition: transport.h:258
void alloc_transport_vectors()
Definition: transport.cc:233
double ** sources_sigma
Definition: transport.h:237
BCField< 3, FieldValue< 3 >::Vector > bc_conc
Definition: transport.h:101
MultiField< 3, FieldValue< 3 >::Scalar > conc_mobile
Calculated concentrations in the mobile zone.
Definition: transport.h:106
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:48
TimeMark::Type target_mark_type
TimeMark type for time marks denoting end of every time interval where transport matrix remains const...
Definition: transport.h:241
void set_initial_condition()
Definition: transport.cc:214
void set_boundary_conditions()
Definition: transport.cc:327
Coupling of a transport model with a reaction model by operator splitting.
void transport_matrix_step_mpi(double time_step)
Definition: mesh.h:108
virtual ~ConvectionTransport()
Definition: transport.cc:173
void get_par_info(int *&el_4_loc, Distribution *&el_ds)
Definition: transport.cc:724
void zero_time_step() override
Definition: transport.cc:468
void output_vector_gather()
Definition: transport.cc:704
Specification of transport model interface.
OutputTime * output_stream_
Definition: transport.h:272
ConvectionTransport(Mesh &init_mesh, const Input::Record &in_rec)
Definition: transport.cc:93
void compute_concentration_sources_for_mass_balance(unsigned int sbi)
Definition: transport.cc:418
double *** out_conc
Definition: transport.h:267
static Input::Type::Selection sorption_type_selection
Definition: transport.h:87
Input::Record output_rec
Record with output specification.
Definition: transport.h:270
void calc_fluxes(vector< vector< double > > &bcd_balance, vector< vector< double > > &bcd_plus_balance, vector< vector< double > > &bcd_minus_balance)
Definition: transport.cc:740
Accessor to the data with type Type::Record.
Definition: accessors.hh:308
The class for outputing data during time.
Definition: output_time.hh:37
Vec * vconc
Concentration vectors for mobile phase.
Definition: transport.h:256
void read_concentration_sources()
void calc_elem_sources(vector< vector< double > > &mass, vector< vector< double > > &src_balance)
Definition: transport.cc:781
FieldSet output_fields
Fields indended for output, i.e. all input fields plus those representing solution.
Definition: transport.h:109
void create_transport_matrix_mpi()
Definition: transport.cc:588
void alloc_transport_structs_mpi()
Definition: transport.cc:279
bool is_convection_matrix_scaled
Definition: transport.h:224
static Input::Type::Selection output_selection
Definition: transport.h:89
VecScatter vconc_out_scatter
Definition: transport.h:247
double ** sources_density
temporary arrays to store constant values of fields over time interval
Definition: transport.h:237
virtual void output_data() override
Write computed fields.
Definition: transport.cc:808
void set_target_time(double target_time)
Definition: transport.cc:567
OutputTime * output_stream()
Definition: transport.h:164
double ** cumulative_corr
Definition: transport.h:264
RegionSet read_boundary_list_item(Input::Record rec)
Override generic method in order to allow specification of the boundary conditions through the old bc...
Class for representation of a vector of fields of the same physical quantity.
Definition: field.hh:309
Distribution * el_ds
Definition: transport.h:277
Template for classes storing finite set of named values.
Field< 3, FieldValue< 3 >::Vector > init_conc
Initial concentrations.
Definition: transport.h:104
void make_transport_partitioning()
Definition: transport.cc:149
double ** sources_conc
Definition: transport.h:237
EqData * get_data()
Definition: transport.h:162