Flow123d  jenkins-Flow123d-windows-release-multijob-285
concentration_model.cc
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 Discontinuous Galerkin method for equation of transport with dispersion.
27  * @author Jan Stebel
28  */
29 
30 #include "input/input_type.hh"
31 #include "mesh/mesh.h"
32 #include "mesh/accessors.hh"
33 #include "flow/darcy_flow_mh.hh"
35 #include "concentration_model.hh"
36 #include "fields/unit_si.hh"
37 
38 
39 
40 using namespace std;
41 using namespace Input::Type;
42 
43 
44 
45 
46 
47 
48 
50 : TransportBase::TransportEqData()
51 {
52  *this+=bc_conc
53  .name("bc_conc")
54  .units( UnitSI().kg().m(-3) )
55  .description("Dirichlet boundary condition (for each substance).")
56  .input_default("0.0")
57  .flags_add( in_rhs );
58  *this+=init_conc
59  .name("init_conc")
60  .units( UnitSI().kg().m(-3) )
61  .description("Initial concentrations.")
62  .input_default("0.0");
63  *this+=disp_l
64  .name("disp_l")
65  .description("Longitudal dispersivity (for each substance).")
66  .units( UnitSI().m() )
67  .input_default("0.0")
69  *this+=disp_t
70  .name("disp_t")
71  .description("Transversal dispersivity (for each substance).")
72  .units( UnitSI().m() )
73  .input_default("0.0")
74  .flags_add( in_main_matrix & in_rhs );
75  *this+=diff_m
76  .name("diff_m")
77  .description("Molecular diffusivity (for each substance).")
78  .units( UnitSI().m(2).s(-1) )
79  .input_default("0.0")
80  .flags_add( in_main_matrix & in_rhs );
81 
82  *this+=output_field
83  .name("conc")
84  .units( UnitSI().kg().m(-3) )
86 }
87 
88 
89 
91 {
92  return data().cross_section.units()*UnitSI().md(1)
93  *data().porosity.units()
94  *data().output_field.units();
95 }
96 
97 
98 IT::Record &ConcentrationTransportModel::get_input_type(const string &implementation, const string &description)
99 {
100  static IT::Record rec = IT::Record(
101  std::string(ModelEqData::name()) + "_" + implementation,
102  description + " for solute transport.")
105  "Names of transported substances.");
106 
107  return rec;
108 }
109 
110 IT::Selection &ConcentrationTransportModel::ModelEqData::get_output_selection_input_type(const string &implementation, const string &description)
111 {
112  static IT::Selection sel = IT::Selection(
113  std::string(ModelEqData::name()) + "_" + implementation + "_Output",
114  "Output record for " + description + " for solute transport.");
115 
116  return sel;
117 }
118 
119 
121  flux_changed(true)
122 {}
123 
124 
126 {
127  substances.initialize(in_rec.val<Input::Array>("substances"));
128 }
129 
130 
132  const ElementAccessor<3> &ele_acc,
133  std::vector<double> &mm_coef)
134 {
135  vector<double> elem_csec(point_list.size()), por_m(point_list.size());
136 
137  data().cross_section.value_list(point_list, ele_acc, elem_csec);
138  data().porosity.value_list(point_list, ele_acc, por_m);
139 
140  for (unsigned int i=0; i<point_list.size(); i++)
141  mm_coef[i] = elem_csec[i]*por_m[i];
142 }
143 
144 
146  double Dm, double alphaL, double alphaT, double porosity, double cross_cut, arma::mat33 &K)
147 {
148  double vnorm = arma::norm(velocity, 2);
149 
150  if (fabs(vnorm) > sqrt(numeric_limits<double>::epsilon()))
151  for (int i=0; i<3; i++)
152  for (int j=0; j<3; j++)
153  K(i,j) = velocity[i]*velocity[j]/(vnorm*vnorm)*(alphaL-alphaT) + alphaT*(i==j?1:0);
154  else
155  K.zeros();
156 
157  // Note that the velocity vector is in fact the Darcian flux,
158  // so to obtain |v| we have to divide vnorm by porosity and cross_section.
159  K = (vnorm*K + (Dm*pow(porosity, 1./3)*porosity*cross_cut)*arma::eye(3,3));
160 }
161 
162 
163 
164 
166  const std::vector<arma::vec3> &velocity,
167  const ElementAccessor<3> &ele_acc,
170 {
171  const unsigned int qsize = point_list.size();
172  const unsigned int n_subst = dif_coef.size();
173  std::vector<arma::vec> Dm(qsize, arma::vec(n_subst) ), alphaL(qsize, arma::vec(n_subst) ), alphaT(qsize, arma::vec(n_subst) );
174  std::vector<double> por_m(qsize), csection(qsize);
175 
176  data().diff_m.value_list(point_list, ele_acc, Dm);
177  data().disp_l.value_list(point_list, ele_acc, alphaL);
178  data().disp_t.value_list(point_list, ele_acc, alphaT);
179  data().porosity.value_list(point_list, ele_acc, por_m);
180  data().cross_section.value_list(point_list, ele_acc, csection);
181 
182  for (unsigned int i=0; i<qsize; i++) {
183  for (unsigned int sbi=0; sbi<n_subst; sbi++) {
184  ad_coef[sbi][i] = velocity[i];
185  calculate_dispersivity_tensor(velocity[i],
186  Dm[i][sbi], alphaL[i][sbi], alphaT[i][sbi], por_m[i], csection[i],
187  dif_coef[sbi][i]);
188  }
189  }
190 }
191 
192 
194  const ElementAccessor<3> &ele_acc,
195  std::vector< arma::vec > &init_values)
196 {
197  data().init_conc.value_list(point_list, ele_acc, init_values);
198 }
199 
200 
202  const ElementAccessor<3> &ele_acc,
203  std::vector< arma::vec > &bc_values)
204 {
205  data().bc_conc.value_list(point_list, ele_acc, bc_values);
206 }
207 
208 
210  const ElementAccessor<3> &ele_acc,
211  std::vector<arma::vec> &sources_value,
212  std::vector<arma::vec> &sources_density,
213  std::vector<arma::vec> &sources_sigma)
214 {
215  const unsigned int qsize = point_list.size();
216  vector<double> csection(qsize);
217  data().cross_section.value_list(point_list, ele_acc, csection);
218  data().sources_conc.value_list(point_list, ele_acc, sources_value);
219  data().sources_density.value_list(point_list, ele_acc, sources_density);
220  data().sources_sigma.value_list(point_list, ele_acc, sources_sigma);
221 
222  for (unsigned int k=0; k<qsize; k++)
223  {
224  sources_density[k] *= csection[k];
225  sources_sigma[k] *= csection[k];
226  }
227 }
228 
229 
231  const ElementAccessor<3> &ele_acc,
232  std::vector<arma::vec> &sources_sigma)
233 {
234  const unsigned int qsize = point_list.size();
235  vector<double> csection(qsize);
236  data().cross_section.value_list(point_list, ele_acc, csection);
237  data().sources_sigma.value_list(point_list, ele_acc, sources_sigma);
238 
239  for (unsigned int k=0; k<qsize; k++)
240  sources_sigma[k] *= csection[k];
241 }
242 
243 
245 {}
246 
247 
248 
Field< 3, FieldValue< 3 >::Vector > disp_l
Longitudal dispersivity (for each substance).
MultiField< 3, FieldValue< 3 >::Scalar > output_field
Accessor to input data conforming to declared Array.
Definition: accessors.hh:558
static Input::Type::AbstractRecord input_type
Common specification of the input record for secondary equations.
Field< 3, FieldValue< 3 >::Vector > disp_t
Transversal dispersivity (for each substance).
static constexpr Mask in_main_matrix
A field is part of main "stiffness matrix" of the equation.
Definition: field_flag.hh:34
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Field< 3, FieldValue< 3 >::Vector > diff_m
Molecular diffusivity (for each substance).
static IT::Record & get_input_type(const string &implementation, const string &description)
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
void compute_dirichlet_bc(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &bc_values) override
static Default obligatory()
Definition: type_record.hh:89
???
void initialize(const Input::Array &in_array)
Read from input array.
Definition: substance.cc:68
Field< 3, FieldValue< 3 >::Vector > sources_conc
Record & derive_from(AbstractRecord &parent)
Definition: type_record.cc:197
Field< 3, FieldValue< 3 >::Scalar > porosity
Mobile porosity.
bool flux_changed
Indicator of change in advection vector field.
static constexpr const char * name()
Field< 3, FieldValue< 3 >::Vector > init_conc
Initial concentrations.
void compute_mass_matrix_coefficient(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &mm_coef) override
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:326
Class for declaration of inputs sequences.
Definition: type_base.hh:239
Specification of transport model interface.
void set_components(SubstanceList &substances, const Input::Record &in_rec) override
Read or set names of solution components.
void compute_init_cond(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &init_values) override
void compute_sources_sigma(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &sources_sigma) override
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:39
static IT::Selection & get_output_selection_input_type(const string &implementation, const string &description)
FieldCommon & input_default(const string &input_default)
Definition: field_common.hh:92
Accessor to the data with type Type::Record.
Definition: accessors.hh:327
void calculate_dispersivity_tensor(const arma::vec3 &velocity, double Dm, double alphaL, double alphaT, double porosity, double cross_cut, arma::mat33 &K)
const Ret val(const string &key) const
BCField< 3, FieldValue< 3 >::Vector > bc_conc
Boundary conditions (Dirichlet) for concentrations.
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh:36
void compute_source_coefficients(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &sources_conc, std::vector< arma::vec > &sources_density, std::vector< arma::vec > &sources_sigma) override
const double epsilon
Definition: mathfce.h:6
Field< 3, FieldValue< 3 >::Vector > sources_sigma
Concentration sources - Robin type, in_flux = sources_sigma * (sources_conc - mobile_conc) ...
FieldCommon & description(const string &description)
Definition: field_common.hh:80
Field< 3, FieldValue< 3 >::Vector > sources_density
Concentration sources - density of substance source, only positive part is used.
FieldCommon & name(const string &name)
Definition: field_common.hh:73
mixed-hybrid model of linear Darcy flow, possibly unsteady.
Record type proxy class.
Definition: type_record.hh:169
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Class for representation SI units of Fields.
Definition: unit_si.hh:31
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 Input::Type::Record input_type
Input type for a substance.
Definition: substance.hh:62
UnitSI & md(int exp=-1)
Method sets value of exponent for m^{-d}, where d is dimension of region.
Definition: unit_si.cc:91
virtual ModelEqData & data()=0
Derived class should implement getter for ModelEqData instance.
Template for classes storing finite set of named values.
Record & declare_key(const string &key, const KeyType &type, const Default &default_value, const string &description)
Definition: type_record.cc:430