Flow123d  JS_before_hm-1-ge159291
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/mesh.h"
21 #include "mesh/accessors.hh"
22 //#include "flow/darcy_flow_mh.hh"
23 
25 #include "concentration_model.hh"
26 #include "tools/unit_si.hh"
27 #include "coupling/balance.hh"
28 
29 
30 
31 using namespace std;
32 using namespace Input::Type;
33 
34 
35 
37  return Selection("Solute_AdvectionDiffusion_BC_Type", "Types of boundary conditions for advection-diffusion solute transport model.")
38  .add_value(bc_inflow, "inflow",
39  "Default transport boundary condition.\n"
40  "On water inflow (($(q_w \\le 0)$)), total flux is given by the reference concentration 'bc_conc'. "
41  "On water outflow we prescribe zero diffusive flux, "
42  "i.e. the mass flows out only due to advection.")
43  .add_value(bc_dirichlet, "dirichlet",
44  "Dirichlet boundary condition (($ c = c_D $)).\n"
45  "The prescribed concentration (($c_D$)) is specified by the field 'bc_conc'.")
46  .add_value(bc_total_flux, "total_flux",
47  "Total mass flux boundary condition.\n"
48  "The prescribed total incoming flux can have the general form (($\\delta(f_N+\\sigma_R(c_R-c) )$)), "
49  "where the absolute flux (($f_N$)) is specified by the field 'bc_flux', "
50  "the transition parameter (($\\sigma_R$)) by 'bc_robin_sigma', "
51  "and the reference concentration (($c_R$)) by 'bc_conc'.")
52  .add_value(bc_diffusive_flux, "diffusive_flux",
53  "Diffusive flux boundary condition.\n"
54  "The prescribed incoming mass flux due to diffusion can have the general form (($\\delta(f_N+\\sigma_R(c_R-c) )$)), "
55  "where the absolute flux (($f_N$)) is specified by the field 'bc_flux', "
56  "the transition parameter (($\\sigma_R$)) by 'bc_robin_sigma', "
57  "and the reference concentration (($c_R$)) by 'bc_conc'.")
58  .close();
59 }
60 
61 
62 
63 
64 
65 
68 {
69  *this+=bc_type
70  .name("bc_type")
71  .description(
72  "Type of boundary condition.")
74  .input_default("\"inflow\"")
77  *this+=bc_dirichlet_value
78  .name("bc_conc")
79  .units( UnitSI().kg().m(-3) )
80  .description("Dirichlet boundary condition (for each substance).")
81  .input_default("0.0")
82  .flags_add( in_rhs );
83  *this+=bc_flux
85  .name("bc_flux")
86  .description("Flux in Neumann boundary condition.")
87  .units( UnitSI().kg().m().s(-1).md() )
88  .input_default("0.0")
89  .flags_add(FieldFlag::in_rhs);
90  *this+=bc_robin_sigma
92  .name("bc_robin_sigma")
93  .description("Conductivity coefficient in Robin boundary condition.")
94  .units( UnitSI().m(4).s(-1).md() )
95  .input_default("0.0")
96  .flags_add(FieldFlag::in_rhs & FieldFlag::in_main_matrix);
97  *this+=init_conc
98  .name("init_conc")
99  .units( UnitSI().kg().m(-3) )
100  .description("Initial values for concentration of substances.")
101  .input_default("0.0");
102  *this+=disp_l
103  .name("disp_l")
104  .description("Longitudinal dispersivity in the liquid (for each substance).")
105  .units( UnitSI().m() )
106  .input_default("0.0")
108  *this+=disp_t
109  .name("disp_t")
110  .description("Transverse dispersivity in the liquid (for each substance).")
111  .units( UnitSI().m() )
112  .input_default("0.0")
113  .flags_add( in_main_matrix & in_rhs );
114  *this+=diff_m
115  .name("diff_m")
116  .description("Molecular diffusivity in the liquid (for each substance).")
117  .units( UnitSI().m(2).s(-1) )
118  .input_default("0.0")
119  .flags_add( in_main_matrix & in_rhs );
120  *this+=rock_density
121  .name("rock_density")
122  .description("Rock matrix density.")
123  .units(UnitSI().kg().m(-3))
124  .input_default("0.0")
126  *this+=sorption_coefficient
127  .name("sorption_coefficient")
128  .description("Coefficient of linear sorption.")
129  .units(UnitSI().m(3).kg(-1))
130  .input_default("0.0")
132 
133  *this+=output_field
134  .name("conc")
135  .description("Concentration solution.")
136  .units( UnitSI().kg().m(-3) )
138 }
139 
140 
141 
142 IT::Record ConcentrationTransportModel::get_input_type(const string &implementation, const string &description)
143 {
144  return IT::Record(
145  std::string(ModelEqData::name()) + "_" + implementation,
146  description + " for solute transport.")
148  .declare_key("solvent_density", IT::Double(0), IT::Default("1.0"),
149  "Density of the solvent [ (($kg.m^{-3}$)) ].");
150 }
151 
153 {
154  // Return empty selection just to provide model specific selection name and description.
155  // The fields are added by TransportDG using an auxiliary selection.
156  return IT::Selection(
157  std::string(ModelEqData::name()) + "_DG_output_fields",
158  "Selection of output fields for Diffusive Solute Transport DG model.");
159 }
160 
161 
163  ConcentrationTransportBase(mesh, in_rec),
164  flux_changed(true),
165  mh_dh(nullptr)
166 {}
167 
168 
170 {
171  solvent_density_ = in_rec.val<double>("solvent_density");
172 }
173 
174 
175 
177  const ElementAccessor<3> &ele_acc,
178  std::vector<double> &mm_coef)
179 {
180  vector<double> elem_csec(point_list.size()), wc(point_list.size()); //por_m(point_list.size()),
181 
182  data().cross_section.value_list(point_list, ele_acc, elem_csec);
183  //data().porosity.value_list(point_list, ele_acc, por_m);
184  data().water_content.value_list(point_list, ele_acc, wc);
185 
186  for (unsigned int i=0; i<point_list.size(); i++)
187  mm_coef[i] = elem_csec[i]*wc[i];
188 }
189 
190 
192  const ElementAccessor<3> &ele_acc,
193  std::vector<std::vector<double> > &ret_coef)
194 {
195  vector<double> elem_csec(point_list.size()),
196  por_m(point_list.size()),
197  rho_l(point_list.size()),
198  rho_s(point_list.size()),
199  sorp_mult(point_list.size());
200 
201  data().cross_section.value_list(point_list, ele_acc, elem_csec);
202  data().porosity.value_list(point_list, ele_acc, por_m);
203  data().rock_density.value_list(point_list, ele_acc, rho_s);
204 
205  // Note: Noe effective water content here, since sorption happen only in the rock (1-porosity).
206  for (unsigned int sbi=0; sbi<substances_.size(); sbi++)
207  {
208  data().sorption_coefficient[sbi].value_list(point_list, ele_acc, sorp_mult);
209  for (unsigned int i=0; i<point_list.size(); i++)
210  {
211  ret_coef[sbi][i] = (1.-por_m[i])*rho_s[i]*sorp_mult[i]*elem_csec[i];
212  }
213  }
214 }
215 
216 
218  const arma::mat33 &Dm, double alphaL, double alphaT, double water_content, double porosity, double cross_cut, arma::mat33 &K)
219 {
220  double vnorm = arma::norm(velocity, 2);
221 
222 
223  // used tortuosity model dues to Millington and Quirk(1961) (should it be with power 10/3 ?)
224  // for an overview of other models see: Chou, Wu, Zeng, Chang (2011)
225  double tortuosity = pow(water_content, 7.0 / 3.0)/ (porosity * porosity);
226 
227  // Note that the velocity vector is in fact the Darcian flux,
228  // so we need not to multiply vnorm by water_content and cross_section.
229  //K = ((alphaL-alphaT) / vnorm) * K + (alphaT*vnorm + Dm*tortuosity*cross_cut*water_content) * arma::eye(3,3);
230 
231  if (fabs(vnorm) > 0) {
232  /*
233  for (int i=0; i<3; i++)
234  for (int j=0; j<3; j++)
235  K(i,j) = (velocity[i]/vnorm)*(velocity[j]);
236  */
237  K = ((alphaL - alphaT) / vnorm) * arma::kron(velocity.t(), velocity);
238 
239  //arma::mat33 abs_diff_mat = arma::abs(K - kk);
240  //double diff = arma::min( arma::min(abs_diff_mat) );
241  //ASSERT( diff < 1e-12 )(diff)(K)(kk);
242  } else
243  K.zeros();
244 
245  // Note that the velocity vector is in fact the Darcian flux,
246  // so to obtain |v| we have to divide vnorm by porosity and cross_section.
247  K += alphaT*vnorm*arma::eye(3,3) + Dm*(tortuosity*cross_cut*water_content);
248 
249 }
250 
251 
252 
253 
255  const std::vector<arma::vec3> &velocity,
256  const ElementAccessor<3> &ele_acc,
259 {
260  const unsigned int qsize = point_list.size();
261  const unsigned int n_subst = dif_coef.size();
262  std::vector<arma::mat33> Dm(qsize);
263  std::vector<double> alphaL(qsize), alphaT(qsize), por_m(qsize), csection(qsize), wc(qsize);
264 
265  data().porosity.value_list(point_list, ele_acc, por_m);
266  data().water_content.value_list(point_list, ele_acc, wc);
267  data().cross_section.value_list(point_list, ele_acc, csection);
268 
269  for (unsigned int sbi=0; sbi<n_subst; sbi++)
270  {
271  data().diff_m[sbi].value_list(point_list, ele_acc, Dm);
272  data().disp_l[sbi].value_list(point_list, ele_acc, alphaL);
273  data().disp_t[sbi].value_list(point_list, ele_acc, alphaT);
274  for (unsigned int i=0; i<qsize; i++)
275  {
276  ad_coef[sbi][i] = velocity[i];
277  calculate_dispersivity_tensor(velocity[i],
278  Dm[i], alphaL[i], alphaT[i], wc[i], por_m[i], csection[i],
279  dif_coef[sbi][i]);
280  }
281  }
282 }
283 
284 
286  const ElementAccessor<3> &ele_acc,
287  std::vector<std::vector<double> > &init_values)
288 {
289  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
290  data().init_conc[sbi].value_list(point_list, ele_acc, init_values[sbi]);
291 }
292 
294  arma::uvec &bc_types)
295 {
296  // Currently the bc types for ConcentrationTransport are numbered in the same way as in TransportDG.
297  // In general we should use some map here.
298  bc_types.resize(n_substances());
299  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
300  bc_types[sbi] = data().bc_type[sbi].value(ele_acc.centre(), ele_acc);
301 }
302 
303 
305  const std::vector<arma::vec3> &point_list,
306  const ElementAccessor<3> &ele_acc,
307  std::vector< double > &bc_flux,
308  std::vector< double > &bc_sigma,
309  std::vector< double > &bc_ref_value)
310 {
311  data().bc_flux[index].value_list(point_list, ele_acc, bc_flux);
312  data().bc_robin_sigma[index].value_list(point_list, ele_acc, bc_sigma);
313  data().bc_dirichlet_value[index].value_list(point_list, ele_acc, bc_ref_value);
314 
315  // Change sign in bc_flux since internally we work with outgoing fluxes.
316  for (auto f : bc_flux) f = -f;
317 }
318 
320  const std::vector<arma::vec3> &point_list,
321  const ElementAccessor<3> &ele_acc,
322  std::vector< double > &bc_sigma)
323 {
324  data().bc_robin_sigma[index].value_list(point_list, ele_acc, bc_sigma);
325 }
326 
327 
329  const ElementAccessor<3> &ele_acc,
330  std::vector<std::vector<double> > &sources_value,
331  std::vector<std::vector<double> > &sources_density,
332  std::vector<std::vector<double> > &sources_sigma)
333 {
334  const unsigned int qsize = point_list.size();
335  vector<double> csection(qsize);
336  data().cross_section.value_list(point_list, ele_acc, csection);
337  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
338  {
339  data().sources_conc[sbi].value_list(point_list, ele_acc, sources_value[sbi]);
340  data().sources_density[sbi].value_list(point_list, ele_acc, sources_density[sbi]);
341  data().sources_sigma[sbi].value_list(point_list, ele_acc, sources_sigma[sbi]);
342 
343  for (unsigned int k=0; k<qsize; k++)
344  {
345  sources_density[sbi][k] *= csection[k];
346  sources_sigma[sbi][k] *= csection[k];
347  }
348  }
349 }
350 
351 
353  const ElementAccessor<3> &ele_acc,
354  std::vector<std::vector<double> > &sources_sigma)
355 {
356  const unsigned int qsize = point_list.size();
357  vector<double> csection(qsize);
358  data().cross_section.value_list(point_list, ele_acc, csection);
359  for (unsigned int sbi=0; sbi<n_substances(); sbi++)
360  {
361  data().sources_sigma[sbi].value_list(point_list, ele_acc, sources_sigma[sbi]);
362 
363  for (unsigned int k=0; k<qsize; k++)
364  sources_sigma[sbi][k] *= csection[k];
365  }
366 }
367 
368 
370 {}
371 
372 
373 void ConcentrationTransportModel::set_balance_object(std::shared_ptr<Balance> balance)
374 {
375  balance_ = balance;
376  subst_idx = balance_->add_quantities(substances_.names());
377 }
378 
379 
380 
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.
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