Flow123d  jenkins-Flow123d-windows32-release-multijob-51
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 
37 
38 
39 using namespace std;
40 using namespace Input::Type;
41 
42 
43 
44 
45 
46 
47 
49 : TransportBase::TransportEqData()
50 {
51  *this+=bc_conc
52  .name("bc_conc")
53  .units( UnitSI().kg().m(-3) )
54  .description("Dirichlet boundary condition (for each substance).")
55  .input_default("0.0")
56  .flags_add( in_rhs );
57  *this+=init_conc
58  .name("init_conc")
59  .units( UnitSI().kg().m(-3) )
60  .description("Initial concentrations.")
61  .input_default("0.0");
62  *this+=disp_l
63  .name("disp_l")
64  .description("Longitudal dispersivity (for each substance).")
65  .units( UnitSI().m() )
66  .input_default("0.0")
68  *this+=disp_t
69  .name("disp_t")
70  .description("Transversal dispersivity (for each substance).")
71  .units( UnitSI().m() )
72  .input_default("0.0")
73  .flags_add( in_main_matrix & in_rhs );
74  *this+=diff_m
75  .name("diff_m")
76  .description("Molecular diffusivity (for each substance).")
77  .units( UnitSI().m(2).s(-1) )
78  .input_default("0.0")
79  .flags_add( in_main_matrix & in_rhs );
80 
81  *this+=output_field
82  .name("conc")
83  .units( UnitSI().kg().m(-3) )
85 }
86 
87 
88 
89 
90 
91 IT::Record &ConcentrationTransportModel::get_input_type(const string &implementation, const string &description)
92 {
93  static IT::Record rec = IT::Record(ModelEqData::name() + "_" + implementation, description + " for solute transport.")
96  "Names of transported substances.");
97 
98  return rec;
99 }
100 
101 IT::Selection &ConcentrationTransportModel::ModelEqData::get_output_selection_input_type(const string &implementation, const string &description)
102 {
103  static IT::Selection sel = IT::Selection(ModelEqData::name() + "_" + implementation + "_Output", "Output record for " + description + " for solute transport.");
104 
105  return sel;
106 }
107 
108 
110  flux_changed(true)
111 {}
112 
113 
115 {
116  in_rec.val<Input::Array>("substances").copy_to(names);
117 }
118 
119 
121  const ElementAccessor<3> &ele_acc,
122  std::vector<double> &mm_coef)
123 {
124  vector<double> elem_csec(point_list.size()), por_m(point_list.size());
125 
126  data().cross_section.value_list(point_list, ele_acc, elem_csec);
127  data().porosity.value_list(point_list, ele_acc, por_m);
128 
129  for (unsigned int i=0; i<point_list.size(); i++)
130  mm_coef[i] = elem_csec[i]*por_m[i];
131 }
132 
133 
135  double Dm, double alphaL, double alphaT, double porosity, double cross_cut, arma::mat33 &K)
136 {
137  double vnorm = arma::norm(velocity, 2);
138 
139  if (fabs(vnorm) > sqrt(numeric_limits<double>::epsilon()))
140  for (int i=0; i<3; i++)
141  for (int j=0; j<3; j++)
142  K(i,j) = velocity[i]*velocity[j]/(vnorm*vnorm)*(alphaL-alphaT) + alphaT*(i==j?1:0);
143  else
144  K.zeros();
145 
146  // Note that the velocity vector is in fact the Darcian flux,
147  // so to obtain |v| we have to divide vnorm by porosity and cross_section.
148  K = (vnorm*K + (Dm*pow(porosity, 1./3)*porosity*cross_cut)*arma::eye(3,3));
149 }
150 
151 
152 
153 
155  const std::vector<arma::vec3> &velocity,
156  const ElementAccessor<3> &ele_acc,
159 {
160  const unsigned int qsize = point_list.size();
161  const unsigned int n_subst = dif_coef.size();
162  std::vector<arma::vec> Dm(qsize, arma::vec(n_subst) ), alphaL(qsize, arma::vec(n_subst) ), alphaT(qsize, arma::vec(n_subst) );
163  std::vector<double> por_m(qsize), csection(qsize);
164 
165  data().diff_m.value_list(point_list, ele_acc, Dm);
166  data().disp_l.value_list(point_list, ele_acc, alphaL);
167  data().disp_t.value_list(point_list, ele_acc, alphaT);
168  data().porosity.value_list(point_list, ele_acc, por_m);
169  data().cross_section.value_list(point_list, ele_acc, csection);
170 
171  for (unsigned int i=0; i<qsize; i++) {
172  for (unsigned int sbi=0; sbi<n_subst; sbi++) {
173  ad_coef[sbi][i] = velocity[i];
174  calculate_dispersivity_tensor(velocity[i],
175  Dm[i][sbi], alphaL[i][sbi], alphaT[i][sbi], por_m[i], csection[i],
176  dif_coef[sbi][i]);
177  }
178  }
179 }
180 
181 
183  const ElementAccessor<3> &ele_acc,
184  std::vector< arma::vec > &init_values)
185 {
186  data().init_conc.value_list(point_list, ele_acc, init_values);
187 }
188 
189 
191  const ElementAccessor<3> &ele_acc,
192  std::vector< arma::vec > &bc_values)
193 {
194  data().bc_conc.value_list(point_list, ele_acc, bc_values);
195 }
196 
197 
199  const ElementAccessor<3> &ele_acc,
200  std::vector<arma::vec> &sources_value,
201  std::vector<arma::vec> &sources_density,
202  std::vector<arma::vec> &sources_sigma)
203 {
204  const unsigned int qsize = point_list.size();
205  vector<double> csection(qsize);
206  data().cross_section.value_list(point_list, ele_acc, csection);
207  data().sources_conc.value_list(point_list, ele_acc, sources_value);
208  data().sources_density.value_list(point_list, ele_acc, sources_density);
209  data().sources_sigma.value_list(point_list, ele_acc, sources_sigma);
210 
211  for (unsigned int k=0; k<qsize; k++)
212  {
213  sources_density[k] *= csection[k];
214  sources_sigma[k] *= csection[k];
215  }
216 }
217 
218 
220  const ElementAccessor<3> &ele_acc,
221  std::vector<arma::vec> &sources_sigma)
222 {
223  const unsigned int qsize = point_list.size();
224  vector<double> csection(qsize);
225  data().cross_section.value_list(point_list, ele_acc, csection);
226  data().sources_sigma.value_list(point_list, ele_acc, sources_sigma);
227 
228  for (unsigned int k=0; k<qsize; k++)
229  sources_sigma[k] *= csection[k];
230 }
231 
232 
234 {}
235 
236 
237 
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:521
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:87
???
Field< 3, FieldValue< 3 >::Vector > sources_conc
Record & derive_from(AbstractRecord &parent)
Definition: type_record.cc:169
Field< 3, FieldValue< 3 >::Scalar > porosity
Mobile porosity.
bool flux_changed
Indicator of change in advection vector field.
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:413
Class for declaration of inputs sequences.
Definition: type_base.hh:230
Specification of transport model interface.
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 non-result fields, that are data fields of an equation.
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:90
Accessor to the data with type Type::Record.
Definition: accessors.hh:308
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:78
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:71
mixed-hybrid model of linear Darcy flow, possibly unsteady.
void set_component_names(std::vector< string > &names, const Input::Record &in_rec) override
Read or set names of solution components.
Record type proxy class.
Definition: type_record.hh:161
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
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:386