Flow123d
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
flow123d
src
transport
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
"
44
#include "
transport/transport_operator_splitting.hh
"
45
46
#include "
fields/field_base.hh
"
47
#include "
fields/field_values.hh
"
48
49
50
class
SorptionImmob
;
51
class
OutputTime
;
52
class
Mesh
;
53
class
Distribution
;
54
class
ConvectionTransport
;
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
*/
82
class
ConvectionTransport
:
public
TransportBase
{
83
public
:
84
85
class
EqData
:
public
TransportBase::TransportEqData
{
86
public
:
87
static
Input::Type::Selection
sorption_type_selection
;
88
89
static
Input::Type::Selection
output_selection
;
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.
95
RegionSet
read_boundary_list_item
(
Input::Record
rec);
96
97
/**
98
* Boundary conditions (Dirichlet) for concentrations.
99
* They are applied only on water inflow part of the boundary.
100
*/
101
BCField<3, FieldValue<3>::Vector
>
bc_conc
;
102
103
/// Initial concentrations.
104
Field<3, FieldValue<3>::Vector
>
init_conc
;
105
//DELETE
106
// Field<3, FieldValue<3>::Scalar> por_imm; ///< Immobile porosity
107
// Field<3, FieldValue<3>::Vector> alpha; ///< Coefficients of non-equilibrium linear mobile-immobile exchange
108
// Field<3, FieldValue<3>::EnumVector> sorp_type; ///< Type of sorption for each substance
109
// Field<3, FieldValue<3>::Vector> sorp_coef0; ///< Coefficient of sorption for each substance
110
// Field<3, FieldValue<3>::Vector> sorp_coef1; ///< Coefficient of sorption for each substance
111
// Field<3, FieldValue<3>::Scalar> phi; ///< solid / solid mobile
112
113
MultiField<3, FieldValue<3>::Scalar
>
conc_mobile
;
///< Calculated concentrations in the mobile zone.
114
115
/// Fields indended for output, i.e. all input fields plus those representing solution.
116
FieldSet
output_fields
;
117
};
118
119
/**
120
* Constructor.
121
*/
122
ConvectionTransport
(
Mesh
&init_mesh,
const
Input::Record
&in_rec);
123
/**
124
* TODO: destructor
125
*/
126
virtual
~ConvectionTransport
();
127
128
/**
129
* Initialize solution at zero time.
130
*/
131
void
zero_time_step
()
override
;
132
/**
133
* Calculates one time step of explicit transport.
134
*/
135
void
update_solution
()
override
;
136
137
/**
138
* Use new flow field vector for construction of convection matrix.
139
* Updates CFL time step constrain.
140
* TODO: Just set the new velocity, postpone update till compute_one_step
141
*/
142
//void set_flow_field_vector(const MH_DofHandler &dh);
143
144
/**
145
* Set cross section of fractures from the Flow equation.
146
*
147
* TODO: Make this and previous part of Transport interface in TransportBase.
148
*/
149
//oid set_cross_section_field(const Field< 3, FieldValue<3>::Scalar > &cross_section);
150
151
152
/**
153
* Set time interval which is considered as one time step by TransportOperatorSplitting.
154
* In particular the velocity field dosn't change over this interval.
155
*
156
* Dependencies:
157
*
158
* velocity, porosity -> matrix, source_vector
159
* matrix -> time_step
160
*
161
* data_read_times -> time_step (not necessary if we won't stick to jump times)
162
* data -> source_vector
163
* time_step -> scaling
164
*
165
*
166
*
167
*/
168
void
set_target_time
(
double
target_time);
169
170
/**
171
* Communicate parallel concentration vectors into sequential output vector.
172
*/
173
void
output_vector_gather
();
//
174
175
/**
176
* @brief Write computed fields.
177
*/
178
virtual
void
output_data
();
179
180
181
/**
182
* Getters.
183
*/
184
inline
EqData
*
get_data
() {
return
&
data_
; }
185
186
inline
OutputTime
*
output_stream
() {
return
output_stream_
; }
187
188
double
***
get_concentration_matrix
();
189
void
get_par_info
(
int
* &
el_4_loc
,
Distribution
* &
el_ds
);
190
int
*
get_el_4_loc
();
191
int
*
get_row_4_el
();
192
193
TimeIntegrationScheme
time_scheme
() {
return
explicit_euler
; }
194
195
private
:
196
197
/**
198
* Assembly convection term part of the matrix and boundary matrix for application of boundary conditions.
199
*
200
* Discretization of the convection term use explicit time scheme and finite volumes with full upwinding.
201
* 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
202
* 3D embedding space)
203
*
204
* In order to get multiplication matrix for explicit transport one have to scale the convection part by the acctual time step and
205
* add time term, i. e. unit matrix (see. transport_matrix_step_mpi)
206
*
207
* Updates CFL time step constrain.
208
*/
209
void
create_transport_matrix_mpi
();
210
211
void
make_transport_partitioning
();
//
212
void
set_initial_condition
();
213
void
read_concentration_sources
();
214
void
set_boundary_conditions
();
215
216
//note: the source of concentration is multiplied by time interval (gives the mass, not the flow like before)
217
void
compute_concentration_sources
(
unsigned
int
sbi);
218
void
compute_concentration_sources_for_mass_balance
(
unsigned
int
sbi);
219
220
/**
221
* Finish explicit transport matrix (time step scaling)
222
*/
223
void
transport_matrix_step_mpi
(
double
time_step);
//
224
225
//DELETE
226
// void transport_dual_porosity( int elm_pos, ElementFullIter elem, int sbi); //
227
// void transport_sorption(int elm_pos, ElementFullIter elem, int sbi); //
228
// void compute_sorption(double conc_avg, double sorp_coef0, double sorp_coef1, unsigned int sorp_type,
229
// double *concx, double *concx_sorb, double Nv, double N); //
230
231
232
void
alloc_transport_vectors
();
233
void
alloc_transport_structs_mpi
();
234
235
/**
236
* Implements the virtual method EquationForMassBalance::calc_fluxes().
237
*/
238
void
calc_fluxes
(
vector
<
vector<double>
> &bcd_balance,
vector
<
vector<double>
> &bcd_plus_balance,
vector
<
vector<double>
> &bcd_minus_balance);
239
/**
240
* Implements the virtual method EquationForMassBalance::calc_elem_sources().
241
*/
242
void
calc_elem_sources
(
vector
<
vector<double>
> &mass,
vector
<
vector<double>
> &src_balance);
243
244
245
/**
246
* Parameters of the equation, some are shared with other implementations since EqData is derived from TransportBase::TransportEqData
247
*/
248
EqData
data_
;
249
250
/**
251
* Indicates if we finished the matrix and add vector by scaling with timestep factor.
252
*/
253
bool
is_convection_matrix_scaled
,
need_time_rescaling
;
254
255
//TODO: remove this and make concentration_matrix only two-dimensional
256
int
sub_problem
;
// 0-only transport,1-transport+dual porosity,
257
// 2-transport+sorption
258
// 3-transport+dual porosity+sorption
259
260
261
double
*
sources_corr
;
262
Vec
v_sources_corr
;
263
264
///temporary arrays to store constant values of fields over time interval
265
//(avoiding calling "field.value()" to often)
266
double
**
sources_density
,
267
**
sources_conc
,
268
**
sources_sigma
;
269
270
TimeMark::Type
target_mark_type
;
///< TimeMark type for time marks denoting end of every time interval where transport matrix remains constant.
271
double
cfl_max_step
;
272
// only local part
273
274
275
276
VecScatter
vconc_out_scatter
;
277
Mat
tm
;
// PETSc transport matrix
278
279
/// Time when the transport matrix was created.
280
/// TODO: when we have our own classes for LA objects, we can use lazy dependence to check
281
/// necessity for matrix update
282
double
transport_matrix_time
;
283
284
/// Concentration vectors for mobile phase.
285
Vec *
vconc
;
// concentration vector
286
/// Concentrations for phase, substance, element
287
double
***
conc
;
288
289
///
290
Vec *
vpconc
;
// previous concentration vector
291
//double ***pconc;
292
Vec *
bcvcorr
;
// boundary condition correction vector
293
Vec *
vcumulative_corr
;
294
double
**
cumulative_corr
;
295
296
Vec *
vconc_out
;
// concentration vector output (gathered)
297
double
***
out_conc
;
298
299
/// Record with output specification.
300
Input::Record
output_rec
;
301
302
OutputTime
*
output_stream_
;
303
304
305
int
*
row_4_el
;
306
int
*
el_4_loc
;
307
Distribution
*
el_ds
;
308
309
friend
class
TransportOperatorSplitting
;
310
};
311
#endif
/* TRANSPORT_H_ */
Generated on Thu May 29 2014 23:14:49 for Flow123d by
1.8.4