Flow123d  JB_transport-112d700
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 
24 #include "concentration_model.hh"
25 #include "tools/unit_si.hh"
26 #include "coupling/balance.hh"
27 #include "fields/field_model.hh"
28 
29 
30 
31 using namespace std;
32 using namespace Input::Type;
33 
34 
35 /*******************************************************************************
36  * Functors of FieldModels
37  */
38 using Sclr = double;
39 using Vect = arma::vec3;
40 using Tens = arma::mat33;
41 
42 // Functor computing velocity norm
43 struct fn_conc_v_norm {
44  inline Sclr operator() (Vect vel) {
45  return arma::norm(vel, 2);
46  }
47 };
48 
49 // Functor computing mass matrix coefficients (cross_section * water_content)
50 struct fn_conc_mass_matrix {
51  inline Sclr operator() (Sclr csec, Sclr wcont) {
52  return csec * wcont;
53  }
54 };
55 
56 // Functor computing retardation coefficients:
57 // (1-porosity) * rock_density * sorption_coefficient * rock_density
58 struct fn_conc_retardation {
59  inline Sclr operator() (Sclr csec, Sclr por_m, Sclr rho_s, Sclr sorp_mult) {
60  return (1.-por_m)*rho_s*sorp_mult*csec;
61  }
62 };
63 
64 // Functor computing sources density output (cross_section * sources_density)
65 struct fn_conc_sources_dens {
66  inline Sclr operator() (Sclr csec, Sclr sdens) {
67  return csec * sdens;
68  }
69 };
70 
71 // Functor computing sources sigma output (cross_section * sources_sigma)
72 struct fn_conc_sources_sigma {
73  inline Sclr operator() (Sclr csec, Sclr ssigma) {
74  return csec * ssigma;
75  }
76 };
77 
78 // Functor computing sources concentration output (sources_conc)
79 struct fn_conc_sources_conc {
80  inline Sclr operator() (Sclr sconc) {
81  return sconc;
82  }
83 };
84 
85 // Functor computing advection coefficient (velocity)
86 struct fn_conc_ad_coef {
87  inline Vect operator() (Vect velocity) {
88  return velocity;
89  }
90 };
91 
92 // Functor computing diffusion coefficient (see notes in function)
94  inline Tens operator() (Tens diff_m, Vect velocity, Sclr v_norm, Sclr alphaL, Sclr alphaT,
95  FMT_UNUSED Sclr water_content, FMT_UNUSED Sclr porosity, Sclr cross_sec) {
96 
97  // used tortuosity model dues to Millington and Quirk(1961) (should it be with power 10/3 ?)
98  // for an overview of other models see: Chou, Wu, Zeng, Chang (2011)
99 
100  // double tortuosity = pow(water_content, 7.0 / 3.0)/ (porosity * porosity);
101 
102  // result
103  Tens K;
104 
105  // Note that the velocity vector is in fact the Darcian flux,
106  // so we need not to multiply vnorm by water_content and cross_section.
107  //K = ((alphaL-alphaT) / vnorm) * K + (alphaT*vnorm + Dm*tortuosity*cross_cut*water_content) * arma::eye(3,3);
108 
109  if (fabs(v_norm) > 0) {
110  /*
111  for (int i=0; i<3; i++)
112  for (int j=0; j<3; j++)
113  K(i,j) = (velocity[i]/vnorm)*(velocity[j]);
114  */
115  K = ((alphaL - alphaT) / v_norm) * arma::kron(velocity.t(), velocity);
116 
117  //arma::mat33 abs_diff_mat = arma::abs(K - kk);
118  //double diff = arma::min( arma::min(abs_diff_mat) );
119  //ASSERT_PERMANENT( diff < 1e-12 )(diff)(K)(kk);
120  } else
121  K.zeros();
122 
123  // Note that the velocity vector is in fact the Darcian flux,
124  // so to obtain |v| we have to divide vnorm by porosity and cross_section.
125  //K += alphaT*v_norm*arma::eye(3,3) + diff_m*(tortuosity*cross_sec*water_content);
126  // We should name it: diffusivity_effective
127  K += alphaT*v_norm*arma::eye(3,3) + diff_m*cross_sec;
128 
129  return K;
130  }
131 };
132 
133 
134 
135 
136 
139 {
140  *this+=bc_type
141  .name("bc_type")
142  .description(
143  "Type of boundary condition.")
145  .input_default("\"inflow\"")
148  *this+=bc_dirichlet_value
149  .name("bc_conc")
150  .units( UnitSI().kg().m(-3) )
151  .description("Dirichlet boundary condition (for each substance).")
152  .input_default("0.0")
153  .flags_add( in_rhs );
154  *this+=bc_flux
156  .name("bc_flux")
157  .description("Flux in Neumann boundary condition.")
158  .units( UnitSI().kg().m().s(-1).md() )
159  .input_default("0.0")
160  .flags_add(FieldFlag::in_rhs);
161  *this+=bc_robin_sigma
163  .name("bc_robin_sigma")
164  .description("Conductivity coefficient in Robin boundary condition.")
165  .units( UnitSI().m(4).s(-1).md() )
166  .input_default("0.0")
168  *this+=init_condition
169  .name("init_conc")
170  .units( UnitSI().kg().m(-3) )
171  .description("Initial values for concentration of substances.")
172  .input_default("0.0");
173  *this+=disp_l
174  .name("disp_l")
175  .description("Longitudinal dispersivity in the liquid (for each substance).")
176  .units( UnitSI().m() )
177  .input_default("0.0")
179  *this+=disp_t
180  .name("disp_t")
181  .description("Transverse dispersivity in the liquid (for each substance).")
182  .units( UnitSI().m() )
183  .input_default("0.0")
185  *this+=diff_m
186  .name("diff_m")
187  .description("Molecular diffusivity in the liquid (for each substance).")
188  .units( UnitSI().m(2).s(-1) )
189  .input_default("0.0")
191  *this+=rock_density
192  .name("rock_density")
193  .description("Rock matrix density.")
194  .units(UnitSI().kg().m(-3))
195  .input_default("0.0")
197  *this+=sorption_coefficient
198  .name("sorption_coefficient")
199  .description("Coefficient of linear sorption.")
200  .units(UnitSI().m(3).kg(-1))
201  .input_default("0.0")
203 
204  *this+=output_field
205  .name("conc")
206  .description("Concentration solution.")
207  .units( UnitSI().kg().m(-3) )
209 
210 
211  // initiaization of FieldModels
212  *this += v_norm.name("v_norm")
213  .description("Velocity norm field.")
214  .input_default("0.0")
215  .units( UnitSI().m().s(-1) );
216 
217  *this += mass_matrix_coef.name("mass_matrix_coef")
218  .description("Matrix coefficients computed by model in mass assemblation.")
219  .input_default("0.0")
220  .units( UnitSI().m(3).md() );
221 
222  *this += retardation_coef.name("retardation_coef")
223  .description("Retardation coefficients computed by model in mass assemblation.")
224  .input_default("0.0")
225  .units( UnitSI().m(3).md() );
226  *this += sources_density_out.name("sources_density_out")
227  .description("Concentration sources output - density of substance source, only positive part is used..")
228  .input_default("0.0")
229  .units( UnitSI().kg().s(-1).md() );
230 
231  *this += sources_sigma_out.name("sources_sigma_out")
232  .description("Concentration sources - Robin type, in_flux = sources_sigma * (sources_conc - mobile_conc).")
233  .input_default("0.0")
234  .units( UnitSI().s(-1).m(3).md() );
235 
236  *this += sources_conc_out.name("sources_conc_out")
237  .description("Concentration sources output.")
238  .input_default("0.0")
239  .units( UnitSI().kg().m(-3) );
240 
241  *this += advection_coef.name("advection_coef")
242  .description("Advection coefficients model.")
243  .input_default("0.0")
244  .units( UnitSI().m().s(-1) );
245 
246  *this += diffusion_coef.name("diffusion_coef")
247  .description("Diffusion coefficients model.")
248  .input_default("0.0")
249  .units( UnitSI().m(2).s(-1) );
250 
251 }
252 
253 
255 {
256  // initialize multifield components
257  sorption_coefficient.setup_components();
258  sources_conc.setup_components();
259  sources_density.setup_components();
260  sources_sigma.setup_components();
261  diff_m.setup_components();
262  disp_l.setup_components();
263  disp_t.setup_components();
264 
265  // create FieldModels
266  v_norm.set(Model<3, FieldValue<3>::Scalar>::create(fn_conc_v_norm(), flow_flux), 0.0);
267  mass_matrix_coef.set(Model<3, FieldValue<3>::Scalar>::create(fn_conc_mass_matrix(), cross_section, water_content), 0.0);
268  retardation_coef.set(
269  Model<3, FieldValue<3>::Scalar>::create_multi(fn_conc_retardation(), cross_section, porosity, rock_density, sorption_coefficient),
270  0.0
271  );
272  sources_density_out.set(Model<3, FieldValue<3>::Scalar>::create_multi(fn_conc_sources_dens(), cross_section, sources_density), 0.0);
273  sources_sigma_out.set(Model<3, FieldValue<3>::Scalar>::create_multi(fn_conc_sources_sigma(), cross_section, sources_sigma), 0.0);
274  sources_conc_out.set(Model<3, FieldValue<3>::Scalar>::create_multi(fn_conc_sources_conc(), sources_conc), 0.0);
275  std::vector<typename Field<3, FieldValue<3>::Vector>::FieldBasePtr> ad_coef_ptr_vec;
276  for (unsigned int sbi=0; sbi<sorption_coefficient.size(); sbi++)
277  ad_coef_ptr_vec.push_back( Model<3, FieldValue<3>::Vector>::create(fn_conc_ad_coef(), flow_flux) );
278  advection_coef.set(ad_coef_ptr_vec, 0.0);
279  diffusion_coef.set(
280  Model<3, FieldValue<3>::Tensor>::create_multi(
281  fn_conc_diff_coef(), diff_m, flow_flux, v_norm, disp_l, disp_t, water_content, porosity, cross_section
282  ),
283  0.0
284  );
285 }
286 
287 
288 
290  return Selection("Solute_AdvectionDiffusion_BC_Type", "Types of boundary conditions for advection-diffusion solute transport model.")
291  .add_value(bc_inflow, "inflow",
292  "Default transport boundary condition.\n"
293  "On water inflow (($(q_w \\le 0)$)), total flux is given by the reference concentration 'bc_conc'. "
294  "On water outflow we prescribe zero diffusive flux, "
295  "i.e. the mass flows out only due to advection.")
296  .add_value(bc_dirichlet, "dirichlet",
297  "Dirichlet boundary condition (($ c = c_D $)).\n"
298  "The prescribed concentration (($c_D$)) is specified by the field 'bc_conc'.")
299  .add_value(bc_total_flux, "total_flux",
300  "Total mass flux boundary condition.\n"
301  "The prescribed total incoming flux can have the general form (($\\delta(f_N+\\sigma_R(c_R-c) )$)), "
302  "where the absolute flux (($f_N$)) is specified by the field 'bc_flux', "
303  "the transition parameter (($\\sigma_R$)) by 'bc_robin_sigma', "
304  "and the reference concentration (($c_R$)) by 'bc_conc'.")
305  .add_value(bc_diffusive_flux, "diffusive_flux",
306  "Diffusive flux boundary condition.\n"
307  "The prescribed incoming mass flux due to diffusion can have the general form (($\\delta(f_N+\\sigma_R(c_R-c) )$)), "
308  "where the absolute flux (($f_N$)) is specified by the field 'bc_flux', "
309  "the transition parameter (($\\sigma_R$)) by 'bc_robin_sigma', "
310  "and the reference concentration (($c_R$)) by 'bc_conc'.")
311  .close();
312 }
313 
314 
315 
316 
318 {
319  // Return empty selection just to provide model specific selection name and description.
320  // The fields are added by TransportDG using an auxiliary selection.
321  return IT::Selection(
322  std::string(ModelEqData::name()) + "_DG_output_fields",
323  "Selection of output fields for Diffusive Solute Transport DG model.");
324 }
325 
326 
327 IT::Record ConcentrationTransportModel::get_input_type(const string &implementation, const string &description)
328 {
329  return IT::Record(
330  std::string(ModelEqData::name()) + "_" + implementation,
331  description + " for solute transport.")
333  .declare_key("solvent_density", IT::Double(0), IT::Default("1.0"),
334  "Density of the solvent [ (($kg.m^{-3}$)) ].");
335 }
336 
337 
340 {}
341 
342 
344 {
345  solvent_density_ = in_rec.val<double>("solvent_density");
346 }
347 
348 
350 {}
351 
352 
353 void ConcentrationTransportModel::set_balance_object(std::shared_ptr<Balance> balance)
354 {
355  balance_ = balance;
356  eq_data().subst_idx_ = balance_->add_quantities(eq_data().substances_.names());
357 }
358 
359 
361 
362 
363 
FieldCommon::units
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
Definition: field_common.hh:153
ConcentrationTransportModel::ModelEqFields::advection_coef
MultiField< 3, FieldValue< 3 >::Vector > advection_coef
Advection coefficients.
Definition: concentration_model.hh:92
fn_conc_ad_coef
Definition: conc_dispersion_model.cc:103
ConcentrationTransportModel::ModelEqData::get_output_selection
static IT::Selection get_output_selection()
Definition: concentration_model.cc:317
FieldFlag::in_rhs
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh:51
ConcentrationTransportModel::ModelEqFields::sources_density_out
MultiField< 3, FieldValue< 3 >::Scalar > sources_density_out
Concentration sources - density output.
Definition: concentration_model.hh:86
ConcentrationTransportModel::ModelEqFields::initialize
void initialize(FMT_UNUSED Input::Record transport_rec)
Definition: concentration_model.cc:254
fn_conc_mass_matrix
Definition: conc_dispersion_model.cc:55
ConcentrationTransportModel::ModelEqFields::init_condition
MultiField< 3, FieldValue< 3 >::Scalar > init_condition
Initial concentrations.
Definition: concentration_model.hh:63
UnitSI::dimensionless
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
ConcentrationTransportModel::set_balance_object
void set_balance_object(std::shared_ptr< Balance > balance) override
Definition: concentration_model.cc:353
ConcentrationTransportModel::ModelEqFields::bc_robin_sigma
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
Transition coefficient in total/diffusive flux b.c.
Definition: concentration_model.hh:61
MultiField::disable_where
auto disable_where(const MultiField< spacedim, typename FieldValue< spacedim >::Enum > &control_field, const vector< FieldEnum > &value_list) -> MultiField &
Definition: multi_field.impl.hh:133
ConcentrationTransportModel::ModelEqData::subst_idx_
vector< unsigned int > subst_idx_
List of indices used to call balance methods for a set of quantities.
Definition: concentration_model.hh:149
ConcentrationTransportModel::eq_data
virtual ModelEqData & eq_data()=0
Derived class should implement getter for ModelEqData instance.
fn_conc_sources_sigma
Definition: conc_dispersion_model.cc:89
Input::Record::val
const Ret val(const string &key) const
Definition: accessors_impl.hh:31
Input::Type::Selection::close
const Selection & close() const
Close the Selection, no more values can be added.
Definition: type_selection.cc:65
ConcentrationTransportModel::ModelEqFields::disp_t
MultiField< 3, FieldValue< 3 >::Scalar > disp_t
Transversal dispersivity (for each substance).
Definition: concentration_model.hh:67
ConcentrationTransportModel::ModelEqFields::diffusion_coef
MultiField< 3, FieldValue< 3 >::Tensor > diffusion_coef
Diffusion coefficients.
Definition: concentration_model.hh:94
Input::Type::Double
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:534
ConcentrationTransportModel::init_from_input
void init_from_input(const Input::Record &in_rec) override
Read necessary data from input record.
Definition: concentration_model.cc:343
fn_conc_sources_conc
Definition: conc_dispersion_model.cc:96
std::vector
Definition: doxy_dummy_defs.hh:7
arma::vec3
Definition: doxy_dummy_defs.hh:17
FieldCommon::flags
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Definition: field_common.hh:192
ConcentrationTransportModel::ModelEqFields::sources_conc_out
MultiField< 3, FieldValue< 3 >::Scalar > sources_conc_out
Concentration sources - concentration output.
Definition: concentration_model.hh:90
FieldCommon::flags_add
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
Definition: field_common.hh:198
ConcentrationTransportModel::ModelEqFields::sources_sigma_out
MultiField< 3, FieldValue< 3 >::Scalar > sources_sigma_out
Concentration sources - sigma output.
Definition: concentration_model.hh:88
EquationBase::mesh
Mesh & mesh()
Definition: equation.hh:181
Input::Type::Default
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Input::Type::Record::derive_from
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
ConcentrationTransportModel::init_balance
void init_balance(const Input::Record &in_rec)
Definition: concentration_model.cc:360
EquationBase::balance
std::shared_ptr< Balance > balance() const
Definition: equation.hh:189
ConcentrationTransportModel::ModelEqFields::get_bc_type_selection
static const Input::Type::Selection & get_bc_type_selection()
Definition: concentration_model.cc:289
accessors.hh
EquationBase::balance_
std::shared_ptr< Balance > balance_
object for calculation and writing the mass balance to file.
Definition: equation.hh:252
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
ConcentrationTransportModel::~ConcentrationTransportModel
~ConcentrationTransportModel() override
Definition: concentration_model.cc:349
FieldFlag::equation_result
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
fn_conc_retardation
Definition: conc_dispersion_model.cc:69
field_model.hh
fn_conc_diff_coef
Definition: concentration_model.cc:93
ConcentrationTransportModel::ModelEqFields::bc_dirichlet
@ bc_dirichlet
Definition: concentration_model.hh:102
ConcentrationTransportModel::ModelEqFields::bc_flux
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_flux
Flux value in total/diffusive flux b.c.
Definition: concentration_model.hh:59
ConcentrationTransportModel::ModelEqFields::bc_type
BCMultiField< 3, FieldValue< 3 >::Enum > bc_type
Type of boundary condition (see also BC_Type)
Definition: concentration_model.hh:55
UnitSI
Class for representation SI units of Fields.
Definition: unit_si.hh:40
ConcentrationTransportModel::ModelEqFields::rock_density
Field< 3, FieldValue< 3 >::Scalar > rock_density
Rock matrix density.
Definition: concentration_model.hh:71
Input::Type::Record::declare_key
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
mesh.h
Input::Type::Selection
Template for classes storing finite set of named values.
Definition: type_selection.hh:65
fn_conc_sources_dens
Definition: conc_dispersion_model.cc:82
ConcentrationTransportModel::ModelEqFields::retardation_coef
MultiField< 3, FieldValue< 3 >::Scalar > retardation_coef
Field represents retardation coefficients due to sorption.
Definition: concentration_model.hh:84
Input::Type
Definition: balance.hh:41
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
Sclr
double Sclr
Definition: field_add_potential.hh:33
FieldFlag::in_main_matrix
static constexpr Mask in_main_matrix
A field is part of main "stiffness matrix" of the equation.
Definition: field_flag.hh:49
ConcentrationTransportModel::ModelEqFields::diff_m
MultiField< 3, FieldValue< 3 >::Tensor > diff_m
Molecular diffusivity (for each substance).
Definition: concentration_model.hh:69
FieldCommon::input_default
FieldCommon & input_default(const string &input_default)
Definition: field_common.hh:140
input_type.hh
ConcentrationTransportModel::ModelEqFields::disp_l
MultiField< 3, FieldValue< 3 >::Scalar > disp_l
Longitudal dispersivity (for each substance).
Definition: concentration_model.hh:65
ConcentrationTransportBase
Definition: transport_operator_splitting.hh:61
Mesh
Definition: mesh.h:362
ConcentrationTransportModel::solvent_density_
double solvent_density_
Density of liquid (a global constant).
Definition: concentration_model.hh:210
Model
Definition: field_model.hh:338
ConcentrationTransportModel::ModelEqFields::bc_inflow
@ bc_inflow
Definition: concentration_model.hh:101
unit_si.hh
ConcentrationTransportModel::ModelEqFields::sorption_coefficient
MultiField< 3, FieldValue< 3 >::Scalar > sorption_coefficient
Coefficient of linear sorption.
Definition: concentration_model.hh:72
ConcentrationTransportModel::ConcentrationTransportModel
ConcentrationTransportModel(Mesh &mesh, const Input::Record &in_rec)
Definition: concentration_model.cc:338
std
Definition: doxy_dummy_defs.hh:5
ConcentrationTransportModel::ModelEqFields::bc_dirichlet_value
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_dirichlet_value
Prescribed concentration for Dirichlet/reference concentration for flux b.c.
Definition: concentration_model.hh:57
ConcentrationTransportModel::ModelEqData::name
static constexpr const char * name()
Definition: concentration_model.hh:123
fn_conc_v_norm
Definition: conc_dispersion_model.cc:48
concentration_model.hh
Discontinuous Galerkin method for equation of transport with dispersion.
ConcentrationTransportBase::get_input_type
static Input::Type::Abstract & get_input_type()
Common specification of the input record for secondary equations.
Definition: transport_operator_splitting.cc:62
balance.hh
ConcentrationTransportModel::ModelEqFields::output_field
MultiField< 3, FieldValue< 3 >::Scalar > output_field
Definition: concentration_model.hh:75
FieldCommon::input_selection
FieldCommon & input_selection(Input::Type::Selection element_selection)
Definition: field_common.hh:172
FieldCommon::description
FieldCommon & description(const string &description)
Definition: field_common.hh:128
TransportEqFields
Definition: transport_operator_splitting.hh:137
Input::Type::Selection::add_value
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.
Definition: type_selection.cc:50
FieldFlag::in_time_term
static constexpr Mask in_time_term
A field is part of time term of the equation.
Definition: field_flag.hh:47
ConcentrationTransportModel::ModelEqFields::ModelEqFields
ModelEqFields()
Definition: concentration_model.cc:137
transport_operator_splitting.hh
ConcentrationTransportModel::ModelEqFields::v_norm
Field< 3, FieldValue< 3 >::Scalar > v_norm
Velocity norm field.
Definition: concentration_model.hh:96
arma::mat33
Definition: doxy_dummy_defs.hh:18
FieldCommon::name
FieldCommon & name(const string &name)
Definition: field_common.hh:121
FieldValue
Definition: field_values.hh:645
FMT_UNUSED
#define FMT_UNUSED
Definition: posix.h:75
ConcentrationTransportModel::ModelEqFields::mass_matrix_coef
Field< 3, FieldValue< 3 >::Scalar > mass_matrix_coef
Field represents coefficients of mass matrix.
Definition: concentration_model.hh:82