Flow123d  release_3.0.0-1191-g7740876
concentration_model.cc
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 concentration_model.cc
15  * @brief Discontinuous Galerkin method for equation of transport with dispersion.
16  * @author Jan Stebel
17  */
18 
19 #include "input/input_type.hh"
20 #include "mesh/side_impl.hh"
21 #include "mesh/mesh.h"
22 #include "mesh/accessors.hh"
23 //#include "flow/darcy_flow_mh.hh"
24 
26 #include "concentration_model.hh"
27 #include "tools/unit_si.hh"
28 #include "coupling/balance.hh"
29 
30 
31 
32 using namespace std;
33 using namespace Input::Type;
34 
35 
36 
38  return Selection("Solute_AdvectionDiffusion_BC_Type", "Types of boundary conditions for advection-diffusion solute transport model.")
39  .add_value(bc_inflow, "inflow",
40  "Default transport boundary condition.\n"
41  "On water inflow (($(q_w \\le 0)$)), total flux is given by the reference concentration 'bc_conc'. "
42  "On water outflow we prescribe zero diffusive flux, "
43  "i.e. the mass flows out only due to advection.")
44  .add_value(bc_dirichlet, "dirichlet",
45  "Dirichlet boundary condition (($ c = c_D $)).\n"
46  "The prescribed concentration (($c_D$)) is specified by the field 'bc_conc'.")
47  .add_value(bc_total_flux, "total_flux",
48  "Total mass flux boundary condition.\n"
49  "The prescribed total incoming flux can have the general form (($\\delta(f_N+\\sigma_R(c_R-c) )$)), "
50  "where the absolute flux (($f_N$)) is specified by the field 'bc_flux', "
51  "the transition parameter (($\\sigma_R$)) by 'bc_robin_sigma', "
52  "and the reference concentration (($c_R$)) by 'bc_conc'.")
53  .add_value(bc_diffusive_flux, "diffusive_flux",
54  "Diffusive flux boundary condition.\n"
55  "The prescribed incoming mass flux due to diffusion can have the general form (($\\delta(f_N+\\sigma_R(c_R-c) )$)), "
56  "where the absolute flux (($f_N$)) is specified by the field 'bc_flux', "
57  "the transition parameter (($\\sigma_R$)) by 'bc_robin_sigma', "
58  "and the reference concentration (($c_R$)) by 'bc_conc'.")
59  .close();
60 }
61 
62 
63 
64 
65 
66 
69 {
70  *this+=bc_type
71  .name("bc_type")
72  .description(
73  "Type of boundary condition.")
75  .input_default("\"inflow\"")
78  *this+=bc_dirichlet_value
79  .name("bc_conc")
80  .units( UnitSI().kg().m(-3) )
81  .description("Dirichlet boundary condition (for each substance).")
82  .input_default("0.0")
83  .flags_add( in_rhs );
84  *this+=bc_flux
86  .name("bc_flux")
87  .description("Flux in Neumann boundary condition.")
88  .units( UnitSI().kg().m().s(-1).md() )
89  .input_default("0.0")
90  .flags_add(FieldFlag::in_rhs);
91  *this+=bc_robin_sigma
93  .name("bc_robin_sigma")
94  .description("Conductivity coefficient in Robin boundary condition.")
95  .units( UnitSI().m(4).s(-1).md() )
96  .input_default("0.0")
97  .flags_add(FieldFlag::in_rhs & FieldFlag::in_main_matrix);
98  *this+=init_conc
99  .name("init_conc")
100  .units( UnitSI().kg().m(-3) )
101  .description("Initial values for concentration of substances.")
102  .input_default("0.0");
103  *this+=disp_l
104  .name("disp_l")
105  .description("Longitudinal dispersivity in the liquid (for each substance).")
106  .units( UnitSI().m() )
107  .input_default("0.0")
109  *this+=disp_t
110  .name("disp_t")
111  .description("Transverse dispersivity in the liquid (for each substance).")
112  .units( UnitSI().m() )
113  .input_default("0.0")
114  .flags_add( in_main_matrix & in_rhs );
115  *this+=diff_m
116  .name("diff_m")
117  .description("Molecular diffusivity in the liquid (for each substance).")
118  .units( UnitSI().m(2).s(-1) )
119  .input_default("0.0")
120  .flags_add( in_main_matrix & in_rhs );
121  *this+=rock_density
122  .name("rock_density")
123  .description("Rock matrix density.")
124  .units(UnitSI().kg().m(-3))
125  .input_default("0.0")
127  *this+=sorption_coefficient
128  .name("sorption_coefficient")
129  .description("Coefficient of linear sorption.")
130  .units(UnitSI().m(3).kg(-1))
131  .input_default("0.0")
133 
134  *this+=output_field
135  .name("conc")
136  .description("Concentration solution.")
137  .units( UnitSI().kg().m(-3) )
139 }
140 
141 
142 
143 IT::Record ConcentrationTransportModel::get_input_type(const string &implementation, const string &description)
144 {
145  return IT::Record(
146  std::string(ModelEqData::name()) + "_" + implementation,
147  description + " for solute transport.")
149  .declare_key("solvent_density", IT::Double(0), IT::Default("1.0"),
150  "Density of the solvent [ (($kg.m^{-3}$)) ].");
151 }
152 
154 {
155  // Return empty selection just to provide model specific selection name and description.
156  // The fields are added by TransportDG using an auxiliary selection.
157  return IT::Selection(
158  std::string(ModelEqData::name()) + "_DG_output_fields",
159  "Selection of output fields for Diffusive Solute Transport DG model.");
160 }
161 
162 
164  ConcentrationTransportBase(mesh, in_rec),
165  flux_changed(true),
166  mh_dh(nullptr)
167 {}
168 
169 
171 {
172  solvent_density_ = in_rec.val<double>("solvent_density");
173 }
174 
175 
176 
178  const ElementAccessor<3> &ele_acc,
179  std::vector<double> &mm_coef)
180 {
181  vector<double> elem_csec(point_list.size()), wc(point_list.size()); //por_m(point_list.size()),
182 
183  data().cross_section.value_list(point_list, ele_acc, elem_csec);
184  //data().porosity.value_list(point_list, ele_acc, por_m);
185  data().water_content.value_list(point_list, ele_acc, wc);
186 
187  for (unsigned int i=0; i<point_list.size(); i++)
188  mm_coef[i] = elem_csec[i]*wc[i];
189 }
190 
191 
193  const ElementAccessor<3> &ele_acc,
194  std::vector<std::vector<double> > &ret_coef)
195 {
196  vector<double> elem_csec(point_list.size()),
197  por_m(point_list.size()),
198  rho_l(point_list.size()),
199  rho_s(point_list.size()),
200  sorp_mult(point_list.size());
201 
202  data().cross_section.value_list(point_list, ele_acc, elem_csec);
203  data().porosity.value_list(point_list, ele_acc, por_m);
204  data().rock_density.value_list(point_list, ele_acc, rho_s);
205 
206  // Note: Noe effective water content here, since sorption happen only in the rock (1-porosity).
207  for (unsigned int sbi=0; sbi<substances_.size(); sbi++)
208  {
209  data().sorption_coefficient[sbi].value_list(point_list, ele_acc, sorp_mult);
210  for (unsigned int i=0; i<point_list.size(); i++)
211  {
212  ret_coef[sbi][i] = (1.-por_m[i])*rho_s[i]*sorp_mult[i]*elem_csec[i];
213  }
214  }
215 }
216 
217 
219  const arma::mat33 &Dm, double alphaL, double alphaT, double water_content, double porosity, double cross_cut, arma::mat33 &K)
220 {
221  double vnorm = arma::norm(velocity, 2);
222 
223 
224  // used tortuosity model dues to Millington and Quirk(1961) (should it be with power 10/3 ?)
225  // for an overview of other models see: Chou, Wu, Zeng, Chang (2011)
226  double tortuosity = pow(water_content, 7.0 / 3.0)/ (porosity * porosity);
227 
228  // Note that the velocity vector is in fact the Darcian flux,
229  // so we need not to multiply vnorm by water_content and cross_section.
230  //K = ((alphaL-alphaT) / vnorm) * K + (alphaT*vnorm + Dm*tortuosity*cross_cut*water_content) * arma::eye(3,3);
231 
232  if (fabs(vnorm) > 0) {
233  /*
234  for (int i=0; i<3; i++)
235  for (int j=0; j<3; j++)
236  K(i,j) = (velocity[i]/vnorm)*(velocity[j]);
237  */
238  K = ((alphaL - alphaT) / vnorm) * arma::kron(velocity.t(), velocity);
239 
240  //arma::mat33 abs_diff_mat = arma::abs(K - kk);
241  //double diff = arma::min( arma::min(abs_diff_mat) );
242  //ASSERT( diff < 1e-12 )(diff)(K)(kk);
243  } else
244  K.zeros();
245 
246  // Note that the velocity vector is in fact the Darcian flux,
247  // so to obtain |v| we have to divide vnorm by porosity and cross_section.
248  K += alphaT*vnorm*arma::eye(3,3) + Dm*(tortuosity*cross_cut*water_content);
249 
250 }
251 
252 
253 
254 
256  const std::vector<arma::vec3> &velocity,
257  const ElementAccessor<3> &ele_acc,
260 {
261  const unsigned int qsize = point_list.size();
262  const unsigned int n_subst = dif_coef.size();
263  std::vector<arma::mat33> Dm(qsize);
264  std::vector<double> alphaL(qsize), alphaT(qsize), por_m(qsize), csection(qsize), wc(qsize);
265 
266  data().porosity.value_list(point_list, ele_acc, por_m);
267  data().water_content.value_list(point_list, ele_acc, wc);
268  data().cross_section.value_list(point_list, ele_acc, csection);
269 
270  for (unsigned int sbi=0; sbi<n_subst; sbi++)
271  {
272  data().diff_m[sbi].value_list(point_list, ele_acc, Dm);
273  data().disp_l[sbi].value_list(point_list, ele_acc, alphaL);
274  data().disp_t[sbi].value_list(point_list, ele_acc, alphaT);
275  for (unsigned int i=0; i<qsize; i++)
276  {
277  ad_coef[sbi][i] = velocity[i];
278  calculate_dispersivity_tensor(velocity[i],
279  Dm[i], alphaL[i], alphaT[i], wc[i], por_m[i], csection[i],
280  dif_coef[sbi][i]);
281  }
282  }
283 }
284 
285 
287  const ElementAccessor<3> &ele_acc,
288  std::vector<std::vector<double> > &init_values)
289 {
290  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
291  data().init_conc[sbi].value_list(point_list, ele_acc, init_values[sbi]);
292 }
293 
295  arma::uvec &bc_types)
296 {
297  // Currently the bc types for ConcentrationTransport are numbered in the same way as in TransportDG.
298  // In general we should use some map here.
299  bc_types.resize(n_substances());
300  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
301  bc_types[sbi] = data().bc_type[sbi].value(ele_acc.centre(), ele_acc);
302 }
303 
304 
306  const std::vector<arma::vec3> &point_list,
307  const ElementAccessor<3> &ele_acc,
308  std::vector< double > &bc_flux,
309  std::vector< double > &bc_sigma,
310  std::vector< double > &bc_ref_value)
311 {
312  data().bc_flux[index].value_list(point_list, ele_acc, bc_flux);
313  data().bc_robin_sigma[index].value_list(point_list, ele_acc, bc_sigma);
314  data().bc_dirichlet_value[index].value_list(point_list, ele_acc, bc_ref_value);
315 
316  // Change sign in bc_flux since internally we work with outgoing fluxes.
317  for (auto f : bc_flux) f = -f;
318 }
319 
321  const std::vector<arma::vec3> &point_list,
322  const ElementAccessor<3> &ele_acc,
323  std::vector< double > &bc_sigma)
324 {
325  data().bc_robin_sigma[index].value_list(point_list, ele_acc, bc_sigma);
326 }
327 
328 
330  const ElementAccessor<3> &ele_acc,
331  std::vector<std::vector<double> > &sources_value,
332  std::vector<std::vector<double> > &sources_density,
333  std::vector<std::vector<double> > &sources_sigma)
334 {
335  const unsigned int qsize = point_list.size();
336  vector<double> csection(qsize);
337  data().cross_section.value_list(point_list, ele_acc, csection);
338  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
339  {
340  data().sources_conc[sbi].value_list(point_list, ele_acc, sources_value[sbi]);
341  data().sources_density[sbi].value_list(point_list, ele_acc, sources_density[sbi]);
342  data().sources_sigma[sbi].value_list(point_list, ele_acc, sources_sigma[sbi]);
343 
344  for (unsigned int k=0; k<qsize; k++)
345  {
346  sources_density[sbi][k] *= csection[k];
347  sources_sigma[sbi][k] *= csection[k];
348  }
349  }
350 }
351 
352 
354  const ElementAccessor<3> &ele_acc,
355  std::vector<std::vector<double> > &sources_sigma)
356 {
357  const unsigned int qsize = point_list.size();
358  vector<double> csection(qsize);
359  data().cross_section.value_list(point_list, ele_acc, csection);
360  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
361  {
362  data().sources_sigma[sbi].value_list(point_list, ele_acc, sources_sigma[sbi]);
363 
364  for (unsigned int k=0; k<qsize; k++)
365  sources_sigma[sbi][k] *= csection[k];
366  }
367 }
368 
369 
371 {}
372 
373 
374 void ConcentrationTransportModel::set_balance_object(std::shared_ptr<Balance> balance)
375 {
376  balance_ = balance;
377  subst_idx = balance_->add_quantities(substances_.names());
378 }
379 
380 
381 
void get_bc_type(const ElementAccessor< 3 > &ele_acc, arma::uvec &bc_types) override
void compute_retardation_coefficient(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &ret_coef) override
SubstanceList substances_
Transported substances.
static const Input::Type::Selection & get_bc_type_selection()
void set_balance_object(std::shared_ptr< Balance > balance) override
ConcentrationTransportModel(Mesh &mesh, const Input::Record &in_rec)
MultiField< 3, FieldValue< 3 >::Scalar > output_field
FieldCommon & input_selection(Input::Type::Selection element_selection)
auto disable_where(const MultiField< spacedim, typename FieldValue< spacedim >::Enum > &control_field, const vector< FieldEnum > &value_list) -> MultiField &
static constexpr Mask in_main_matrix
A field is part of main "stiffness matrix" of the equation.
Definition: field_flag.hh:49
MultiField< 3, FieldValue< 3 >::Scalar > init_conc
Initial concentrations.
const std::vector< std::string > & names()
Definition: substance.hh:85
void init_from_input(const Input::Record &in_rec) override
Read necessary data from input record.
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
MultiField< 3, FieldValue< 3 >::Scalar > sorption_coefficient
Coefficient of linear sorption.
void calculate_dispersivity_tensor(const arma::vec3 &velocity, const arma::mat33 &Dm, double alphaL, double alphaT, double water_content, double porosity, double cross_cut, arma::mat33 &K)
void get_flux_bc_data(unsigned int index, const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &bc_flux, std::vector< double > &bc_sigma, std::vector< double > &bc_ref_value) override
Return data for diffusive or total flux b.c.
MultiField< 3, FieldValue< 3 >::Scalar > sources_density
Concentration sources - density of substance source, only positive part is used.
void compute_init_cond(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &init_values) override
Definition: mesh.h:76
MultiField< 3, FieldValue< 3 >::Scalar > sources_conc
bool flux_changed
Indicator of change in advection vector field.
static constexpr const char * name()
void compute_mass_matrix_coefficient(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &mm_coef) override
static Input::Type::Abstract & get_input_type()
Common specification of the input record for secondary equations.
Discontinuous Galerkin method for equation of transport with dispersion.
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
virtual void value_list(const std::vector< Point > &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list) const
Definition: field.hh:403
unsigned int size() const
Definition: substance.hh:87
double solvent_density_
Density of liquid (a global constant).
Field< 3, FieldValue< 3 >::Scalar > rock_density
Rock matrix density.
MultiField< 3, FieldValue< 3 >::Scalar > disp_t
Transversal dispersivity (for each substance).
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
Definition: accessors.hh:285
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
BCMultiField< 3, FieldValue< 3 >::Enum > bc_type
Type of boundary condition (see also BC_Type)
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_dirichlet_value
Prescribed concentration for Dirichlet/reference concentration for flux b.c.
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:541
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
Mesh & mesh()
Definition: equation.hh:176
static constexpr Mask in_time_term
A field is part of time term of the equation.
Definition: field_flag.hh:47
FieldCommon & input_default(const string &input_default)
MultiField< 3, FieldValue< 3 >::Scalar > sources_sigma
Concentration sources - Robin type, in_flux = sources_sigma * (sources_conc - mobile_conc) ...
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
const Ret val(const string &key) const
Selection & add_value(const int value, const std::string &key, const std::string &description="", TypeBase::attribute_map attributes=TypeBase::attribute_map())
Adds one new value with name given by key to the Selection.
Field< 3, FieldValue< 3 >::Scalar > porosity
Mobile porosity - usually saturated water content in the case of unsaturated flow model...
std::shared_ptr< Balance > balance_
object for calculation and writing the mass balance to file.
Definition: equation.hh:235
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:503
unsigned int n_substances() override
Returns number of transported substances.
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh:51
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
Transition coefficient in total/diffusive flux b.c.
Field< 3, FieldValue< 3 >::Scalar > water_content
Water content - result of unsaturated water flow model or porosity.
FieldCommon & description(const string &description)
MultiField< 3, FieldValue< 3 >::Scalar > disp_l
Longitudal dispersivity (for each substance).
void get_flux_bc_sigma(unsigned int index, const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &bc_sigma) override
Return transition coefficient for flux b.c.
const Selection & close() const
Close the Selection, no more values can be added.
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_flux
Flux value in total/diffusive flux b.c.
FieldCommon & name(const string &name)
MultiField< 3, FieldValue< 3 >::TensorFixed > diff_m
Molecular diffusivity (for each substance).
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
Record type proxy class.
Definition: type_record.hh:182
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Class for representation SI units of Fields.
Definition: unit_si.hh:40
void compute_advection_diffusion_coefficients(const std::vector< arma::vec3 > &point_list, const std::vector< arma::vec3 > &velocity, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< arma::vec3 > > &ad_coef, std::vector< std::vector< arma::mat33 > > &dif_coef) override
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
virtual ModelEqData & data()=0
Derived class should implement getter for ModelEqData instance.
void compute_sources_sigma(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &sources_sigma) override
Template for classes storing finite set of named values.
void compute_source_coefficients(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &sources_conc, std::vector< std::vector< double > > &sources_density, std::vector< std::vector< double > > &sources_sigma) override