Flow123d  JS_before_hm-984-g3a19f2f
sorption_base.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 sorption_base.cc
15  * @brief
16  */
17 
23 #include "reaction/isotherm.hh"
24 
25 #include "system/system.hh"
26 #include "system/sys_profiler.hh"
27 
28 #include "la/distribution.hh"
29 #include "mesh/mesh.h"
30 #include "mesh/elements.h"
31 #include "mesh/accessors.hh"
32 #include "input/type_selection.hh"
33 
34 #include "fields/field_set.hh"
35 #include "fields/field_fe.hh"
37 
38 using namespace Input::Type;
39 
41  return Selection("SorptionType")
42  .add_value(Isotherm::none,"none", "No sorption considered.")
43  .add_value(Isotherm::linear, "linear",
44  "Linear isotherm runs the concentration exchange between liquid and solid.")
45  .add_value(Isotherm::langmuir, "langmuir",
46  "Langmuir isotherm runs the concentration exchange between liquid and solid.")
47  .add_value(Isotherm::freundlich, "freundlich",
48  "Freundlich isotherm runs the concentration exchange between liquid and solid.")
49  .close();
50 }
51 
52 
53 
55  return Record("SorptionAux", "AUXILIARY RECORD. Should not be directly part of the input tree.")
56  .declare_key("substances", Array(String(),1), Default::obligatory(),
57  "Names of the substances that take part in the sorption model.")
58  .declare_key("solvent_density", Double(0.0), Default("1.0"),
59  "Density of the solvent.")
60  .declare_key("substeps", Integer(1), Default("1000"),
61  "Number of equidistant substeps, molar mass and isotherm intersections")
62  .declare_key("solubility", Array(Double(0.0)), Default::optional(), //("-1.0"), //
63  "Specifies solubility limits of all the sorbing species.")
64  .declare_key("table_limits", Array(Double(-1.0)), Default::optional(), //("-1.0"), //
65  "Specifies the highest aqueous concentration in the isotherm function interpolation table. "
66  "Use any negative value for an automatic choice according to current maximal concentration (default and recommended). "
67  "Use '0' to always evaluate isotherm function directly (can be very slow). "
68  "Use a positive value to set the interpolation table limit manually "
69  "(if aqueous concentration is higher, then the isotherm function is evaluated directly).")
70  .declare_key("input_fields", Array(EqData("","").input_data_set_.make_field_descriptor_type("Sorption")), Default::obligatory(), //
71  "Containes region specific data necessary to construct isotherms.")//;
72  .declare_key("reaction_liquid", ReactionTerm::it_abstract_reaction(), Default::optional(), "Reaction model following the sorption in the liquid.")
73  .declare_key("reaction_solid", ReactionTerm::it_abstract_reaction(), Default::optional(), "Reaction model following the sorption in the solid.")
74  .close();
75 }
76 
77 
78 SorptionBase::EqData::EqData(const string &output_field_name, const string &output_field_desc)
79 {
80  *this += rock_density.name("rock_density")
81  .description("Rock matrix density.")
82  .input_default("0.0")
83  .units( UnitSI().kg().m(-3) );
84 
85  *this += sorption_type.name("sorption_type")
86  .description("Considered sorption is described by selected isotherm.\n"
87  "If porosity on an element is equal to 1.0 (or even higher), meaning no sorbing surface, then type 'none' will be selected automatically.")
88  .input_selection(get_sorption_type_selection())
89  .units( UnitSI::dimensionless() );
90 
91  *this += distribution_coefficient.name("distribution_coefficient")
92  .description("Distribution coefficient (( $k_l, k_F, k_L $)) of linear, Freundlich or Langmuir isotherm respectively.")
93  .input_default("1.0")
94  .units( UnitSI().m(3).kg(-1) );
95 
96  *this += isotherm_other.name("isotherm_other")
97  .description("Additional parameter (($ \\alpha $)) of nonlinear isotherms.")
98  .input_default("1.0")
99  .units( UnitSI::dimensionless() );
100 
101  *this += init_conc_solid.name("init_conc_solid")
102  .description("Initial solid concentration of substances. It is a vector: one value for every substance.")
103  .input_default("0")
104  .units( UnitSI().dimensionless() );
105 
106  input_data_set_ += *this;
107 
108  // porosity field is set from governing equation (transport) later
109  // hence we do not add it to the input_data_set_
110  *this += porosity
111  .name("porosity")
112  .units( UnitSI::dimensionless() )
113  .flags(FieldFlag::input_copy)
114  .set_limits(0.0);
115 
116  output_fields += *this;
117  output_fields += conc_solid.name(output_field_name)
118  .description(output_field_desc)
119  .units( UnitSI().dimensionless() )
121 }
122 
123 
125  : ReactionTerm(init_mesh, in_rec),
126  data_(nullptr)
127 {
128  // creating reaction from input and setting their parameters
129  make_reactions();
130 }
131 
132 
134 {
135  if (data_ != nullptr) delete data_;
136 }
137 
139 {
141 
142  reactions_it = input_record_.find<Input::AbstractRecord>("reaction_liquid");
143  if ( reactions_it )
144  {
145  // TODO: allowed instances in this case are only
146  // FirstOrderReaction, RadioactiveDecay
147  reaction_liquid = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
148  } else
149  {
150  reaction_liquid = nullptr;
151  }
152 
153  reactions_it = input_record_.find<Input::AbstractRecord>("reaction_solid");
154  if ( reactions_it )
155  {
156  // TODO: allowed instances in this case are only
157  // FirstOrderReaction, RadioactiveDecay
158  reaction_solid = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
159  } else
160  {
161  reaction_solid = nullptr;
162  }
163 }
164 
166 {
167  ASSERT_PTR(time_).error("Time governor has not been set yet.\n");
168  ASSERT(output_stream_).error("Null output stream.\n");
169  ASSERT_LT(0, substances_.size());
170 
171  initialize_substance_ids(); //computes present substances and sets indices
172  initialize_from_input(); //reads non-field data from input
173 
174  //isotherms array resized bellow
175  unsigned int nr_of_regions = mesh_->region_db().bulk_size();
176  isotherms.resize(nr_of_regions);
177  max_conc.resize(nr_of_regions);
178  for(unsigned int i_reg = 0; i_reg < nr_of_regions; i_reg++)
179  {
180  isotherms[i_reg].resize(n_substances_);
181  max_conc[i_reg].resize(n_substances_, 0.0);
182  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++)
183  {
184  isotherms[i_reg][i_subst] = Isotherm();
185  }
186  }
187 
189 
190  if(reaction_liquid)
191  {
192  reaction_liquid->substances(substances_)
193  .concentration_fields(conc_mobile_fe)
194  .set_time_governor(*time_);
195  reaction_liquid->initialize();
196  }
197  if(reaction_solid)
198  {
199  reaction_solid->substances(substances_)
200  .concentration_fields(data_->conc_solid_fe)
201  .set_time_governor(*time_);
202  reaction_solid->initialize();
203  }
204 }
205 
206 
208 {
209  Input::Array substances_array = input_record_.val<Input::Array>("substances");
210  unsigned int k, global_idx, i_subst = 0;
211  bool found;
212  Input::Iterator<string> spec_iter = substances_array.begin<string>();
213 
214  for(; spec_iter != substances_array.end(); ++spec_iter, i_subst++)
215  {
216  //finding the name of a substance in the global array of names
217  found = false;
218  for(k = 0; k < substances_.size(); k++)
219  {
220  if (*spec_iter == substances_[k].name())
221  {
222  global_idx = k;
223  found = true;
224  break;
225  }
226  }
227 
228  if(!found)
229  THROW(ReactionTerm::ExcUnknownSubstance()
230  << ReactionTerm::EI_Substance(*spec_iter)
231  << substances_array.ei_address());
232 
233  //finding the global index of substance in the local array
234  found = false;
235  for(k = 0; k < substance_global_idx_.size(); k++)
236  {
237  if(substance_global_idx_[k] == global_idx)
238  {
239  found = true;
240  break;
241  }
242  }
243 
244  if(!found)
245  {
246  substance_global_idx_.push_back(global_idx);
247  }
248 
249  }
251 }
252 
254 {
255  // read number of interpolation steps - value checked by the record definition
256  n_interpolation_steps_ = input_record_.val<int>("substeps");
257 
258  // read the density of solvent - value checked by the record definition
259  solvent_density_ = input_record_.val<double>("solvent_density");
260 
261  // read the solubility vector
262  Input::Iterator<Input::Array> solub_iter = input_record_.find<Input::Array>("solubility");
263  if( solub_iter )
264  {
265  if (solub_iter->Array::size() != n_substances_)
266  {
267  THROW(SorptionBase::ExcSubstanceCountMatch()
268  << SorptionBase::EI_ArrayName("solubility")
270  // there is no way to get ei_address from 'solub_iter', only as a string
271  }
272 
273  else solub_iter->copy_to(solubility_vec_);
274  }
275  else{
276  // fill solubility_vec_ with zeros
277  solubility_vec_.clear();
278  solubility_vec_.resize(n_substances_,0.0);
279  }
280 
281  // read the interpolation table limits
282  Input::Iterator<Input::Array> interp_table_limits = input_record_.find<Input::Array>("table_limits");
283  if( interp_table_limits )
284  {
285  if (interp_table_limits->Array::size() != n_substances_)
286  {
287  THROW(SorptionBase::ExcSubstanceCountMatch()
288  << SorptionBase::EI_ArrayName("table_limits")
290  // there is no way to get ei_address from 'interp_table_limits', only as a string
291  }
292 
293  else interp_table_limits->copy_to(table_limit_);
294  }
295  else{
296  // fill table_limit_ with negative values -> automatic limit
297  table_limit_.clear();
298  table_limit_.resize(n_substances_,-1.0);
299  }
300 }
301 
303 {
304  ASSERT_GT(n_substances_, 0).error("Number of substances is wrong, they might have not been set yet.\n");
305 
306  // create vector of substances that are involved in sorption
307  // and initialize data_ with their names
308  std::vector<std::string> substances_sorption;
309  for (unsigned int i : substance_global_idx_)
310  substances_sorption.push_back(substances_[i].name());
311  data_->set_components(substances_sorption);
312 
313  // read fields from input file
315 
316  data_->set_mesh(*mesh_);
317 
318  //initialization of output
319  //output_array = input_record_.val<Input::Array>("output_fields");
324 
325  //creating field fe and output multifield for sorbed concentrations
326  data_->conc_solid_fe.resize(substances_.size());
327  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
328  {
329  data_->conc_solid_fe[sbi] = create_field_fe< 3, FieldValue<3>::Scalar >(dof_handler_);
330  data_->conc_solid[sbi].set_field(mesh_->region_db().get_region_set("ALL"), data_->conc_solid_fe[sbi], 0);
331  }
332  //output_stream_->add_admissible_field_names(output_array);
334 }
335 
336 
338 {
339  ASSERT_PTR(time_).error("Time governor has not been set yet.\n");
340  ASSERT(output_stream_).error("Null output stream.\n");
341  ASSERT_LT(0, substances_.size());
342 
344  std::stringstream ss; // print warning message with table of uninitialized fields
345  if ( FieldCommon::print_message_table(ss, "sorption") ) {
346  WarningOut() << ss.str();
347  }
349 
350  update_max_conc();
351  make_tables();
352 
353  // write initial condition
354  //data_->output_fields.set_time(time_->step(), LimitSide::right);
355  //data_->output_fields.output(output_stream_);
356 
357  if(reaction_liquid) reaction_liquid->zero_time_step();
358  if(reaction_solid) reaction_solid->zero_time_step();
359 
360  output_data();
361 }
362 
364 {
365  for ( DHCellAccessor dh_cell : dof_handler_->own_range() ) {
366  IntIdx dof_p0 = dh_cell.get_loc_dof_indices()[0];
367  const ElementAccessor<3> ele = dh_cell.elm();
368 
369  //setting initial solid concentration for substances involved in adsorption
370  for (unsigned int sbi = 0; sbi < n_substances_; sbi++)
371  {
372  int subst_id = substance_global_idx_[sbi];
373  data_->conc_solid_fe[subst_id]->vec()[dof_p0] = data_->init_conc_solid[sbi].value(ele.centre(), ele);
374  }
375  }
376 }
377 
378 
380 {
381  data_->set_time(time_->step(), LimitSide::right); // set to the last computed time
382 
383  // if parameters changed during last time step, reinit isotherms and eventualy
384  // update interpolation tables in the case of constant rock matrix parameters
385  make_tables();
386  clear_max_conc();
387 
388  START_TIMER("Sorption");
389  for ( DHCellAccessor dh_cell : dof_handler_->own_range() )
390  {
391  compute_reaction(dh_cell);
392  }
393  END_TIMER("Sorption");
394 
395  if(reaction_liquid) reaction_liquid->update_solution();
396  if(reaction_solid) reaction_solid->update_solution();
397 }
398 
399 void SorptionBase::isotherm_reinit(unsigned int i_subst, const ElementAccessor<3> &elem)
400 {
401  START_TIMER("SorptionBase::isotherm_reinit");
402 
403  double mult_coef = data_->distribution_coefficient[i_subst].value(elem.centre(),elem);
404  double second_coef = data_->isotherm_other[i_subst].value(elem.centre(),elem);
405 
406  int reg_idx = elem.region().bulk_idx();
407  Isotherm & isotherm = isotherms[reg_idx][i_subst];
408 
409  bool limited_solubility_on = solubility_vec_[i_subst] > 0.0;
410 
411  // in case of no sorbing surface, set isotherm type None
412  if( common_ele_data.no_sorbing_surface_cond <= std::numeric_limits<double>::epsilon())
413  {
414  isotherm.reinit(Isotherm::none, false, solvent_density_,
416  0,0,0);
417  return;
418  }
419 
420  if ( common_ele_data.scale_sorbed <= 0.0)
421  xprintf(UsrErr, "Scaling parameter in sorption is not positive. Check the input for rock density and molar mass of %d. substance.",i_subst);
422 
423  isotherm.reinit(Isotherm::SorptionType(data_->sorption_type[i_subst].value(elem.centre(),elem)),
424  limited_solubility_on, solvent_density_,
426  solubility_vec_[i_subst], mult_coef, second_coef);
427 }
428 
430 {
431  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++)
432  {
433  isotherm_reinit(i_subst, elem);
434  }
435 }
436 
438 {
439  unsigned int reg_idx, i_subst;
440 
441  // clear max concetrations array
442  unsigned int nr_of_regions = mesh_->region_db().bulk_size();
443  for(reg_idx = 0; reg_idx < nr_of_regions; reg_idx++)
444  for(i_subst = 0; i_subst < n_substances_; i_subst++)
445  max_conc[reg_idx][i_subst] = 0.0;
446 }
447 
449 {
450  unsigned int reg_idx, i_subst, subst_id;
451 
452  clear_max_conc();
453 
454  for ( DHCellAccessor dh_cell : dof_handler_->own_range() ) {
455  IntIdx dof_p0 = dh_cell.get_loc_dof_indices()[0];
456  reg_idx = dh_cell.elm().region().bulk_idx();
457  for(i_subst = 0; i_subst < n_substances_; i_subst++){
458  subst_id = substance_global_idx_[i_subst];
459  max_conc[reg_idx][i_subst] = std::max(max_conc[reg_idx][i_subst], conc_mobile_fe[subst_id]->vec()[dof_p0]);
460  }
461  }
462 }
463 
465 {
466  START_TIMER("SorptionBase::make_tables");
467  try
468  {
469  ElementAccessor<3> elm;
470  for(const Region &reg_iter: this->mesh_->region_db().get_region_set("BULK"))
471  {
472  int reg_idx = reg_iter.bulk_idx();
473  // true if data has been changed and are constant on the region
474  bool call_reinit = data_->changed() && data_->is_constant(reg_iter);
475 
476  if(call_reinit)
477  {
478  ElementAccessor<3> elm(this->mesh_, reg_iter);
479 // DebugOut().fmt("isotherm reinit\n");
481  isotherm_reinit_all(elm);
482  }
483 
484  // find table limit and create interpolation table for every substance
485  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++){
486 
487  // clear interpolation tables, if not spacially constant OR switched off
488  if(! data_->is_constant(reg_iter) || table_limit_[i_subst] == 0.0){
489  isotherms[reg_idx][i_subst].clear_table();
490 // DebugOut().fmt("limit: 0.0 -> clear table\n");
491  continue;
492  }
493 
494  // if true then make_table will be called at the end
495  bool call_make_table = call_reinit;
496  // initialy try to keep the current table limit (it is zero at zero time step)
497  double subst_table_limit = isotherms[reg_idx][i_subst].table_limit();
498 
499  // if automatic, possibly remake tables with doubled range when table maximum was reached
500  if(table_limit_[i_subst] < 0.0)
501  {
502  if(subst_table_limit < max_conc[reg_idx][i_subst])
503  {
504  call_make_table = true;
505  subst_table_limit = 2*max_conc[reg_idx][i_subst];
506 // DebugOut().fmt("limit: max conc\n");
507  }
508  }
509  // if not automatic, set given table limit
510  else
511  {
512  subst_table_limit = table_limit_[i_subst];
513  }
514 
515  if(call_make_table){
516  isotherms[reg_idx][i_subst].make_table(n_interpolation_steps_, subst_table_limit);
517 // DebugOut().fmt("reg: {} i_subst {}: table_limit = {}\n", reg_idx, i_subst, isotherms[reg_idx][i_subst].table_limit());
518  }
519  }
520  }
521  }
522  catch(ExceptionBase const &e)
523  {
524  e << input_record_.ei_address();
525  throw;
526  }
527 }
528 
530 {
531  const ElementAccessor<3> ele = dh_cell.elm();
532  int reg_idx = ele.region().bulk_idx();
533  IntIdx dof_p0 = dh_cell.get_loc_dof_indices()[0];
534  unsigned int i_subst, subst_id;
535  // for checking, that the common element data are computed once at maximum
536  bool is_common_ele_data_valid = false;
537 
538  try{
539  for(i_subst = 0; i_subst < n_substances_; i_subst++)
540  {
541  subst_id = substance_global_idx_[i_subst];
542  Isotherm & isotherm = isotherms[reg_idx][i_subst];
543  if (isotherm.is_precomputed()){
544 // DebugOut().fmt("isotherms precomputed - interpolate, subst[{}]\n", i_subst);
545  isotherm.interpolate(conc_mobile_fe[subst_id]->vec()[dof_p0],
546  data_->conc_solid_fe[subst_id]->vec()[dof_p0]);
547  }
548  else{
549 // DebugOut().fmt("isotherms reinit - compute , subst[{}]\n", i_subst);
550  if(! is_common_ele_data_valid){
552  is_common_ele_data_valid = true;
553  }
554 
555  isotherm_reinit(i_subst, ele);
556  isotherm.compute(conc_mobile_fe[subst_id]->vec()[dof_p0],
557  data_->conc_solid_fe[subst_id]->vec()[dof_p0]);
558  }
559 
560  // update maximal concentration per region (optimization for interpolation)
561  if(table_limit_[i_subst] < 0)
562  max_conc[reg_idx][i_subst] = std::max(max_conc[reg_idx][i_subst],
563  conc_mobile_fe[subst_id]->vec()[dof_p0]);
564  }
565  }
566  catch(ExceptionBase const &e)
567  {
568  e << input_record_.ei_address();
569  throw;
570  }
571 }
572 
573 
574 /**************************************** OUTPUT ***************************************************/
575 
577 {
579  // Register fresh output data
580  data_->output_fields.output(time().step());
581 }
TimeGovernor & time()
Definition: equation.hh:149
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:214
std::shared_ptr< ReactionTerm > reaction_solid
Iterator< ValueType > begin() const
MultiField< 3, FieldValue< 3 >::Scalar > conc_solid
Calculated sorbed concentrations, for output only.
RegionSet get_region_set(const std::string &set_name) const
Definition: region.cc:328
void make_reactions()
Accessor to input data conforming to declared Array.
Definition: accessors.hh:566
double solvent_density_
ArmaVec< double, N > vec
Definition: armor.hh:861
EI_Address ei_address() const
Definition: accessors.cc:314
void compute(double &c_aqua, double &c_sorbed)
Definition: isotherm.cc:68
std::vector< std::vector< Isotherm > > isotherms
virtual ~SorptionBase(void)
void initialize_substance_ids()
Reads names of substances from input and creates indexing to global vector of substance.
unsigned int bulk_size() const
Definition: region.cc:268
const std::vector< std::string > & names()
Definition: substance.hh:85
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
int IntIdx
Definition: index_types.hh:25
void output(TimeStep step)
#define ASSERT_GT(a, b)
Definition of comparative assert macro (Greater Than)
Definition: asserts.hh:312
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
SubstanceList substances_
Names belonging to substances.
void zero_time_step() override
Definition: mesh.h:78
Iterator< Ret > find(const string &key) const
Cell accessor allow iterate over DOF handler cells.
EquationOutput output_fields
Fields indended for output, i.e. all input fields plus those representing solution.
static const Input::Type::Selection & get_sorption_type_selection()
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
void compute_reaction(const DHCellAccessor &dh_cell) override
Compute reaction on a single element.
const RegionDB & region_db() const
Definition: mesh.h:143
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
const TimeStep & step(int index=-1) const
std::vector< double > table_limit_
Class for declaration of the integral input data.
Definition: type_base.hh:483
static Input::Type::Abstract & it_abstract_reaction()
void reinit(enum SorptionType sorption_type, bool limited_solubility_on, double aqua_density, double scale_aqua, double scale_sorbed, double c_aqua_limit, double mult_coef, double second_coef)
Definition: isotherm.cc:44
unsigned int size() const
Definition: substance.hh:87
Class for declaration of inputs sequences.
Definition: type_base.hh:339
std::shared_ptr< OutputTime > output_stream_
Pointer to a transport output stream.
ReactionTerm(Mesh &init_mesh, Input::Record in_rec)
Constructor.
Class ReactionTerm is an abstract class representing reaction term in transport.
std::vector< unsigned int > substance_global_idx_
Mapping from local indexing of substances to global.
void initialize_fields()
Initializes field sets.
void isotherm_reinit(unsigned int i_subst, const ElementAccessor< 3 > &elm)
Reinitializes the isotherm.
Class Dual_por_exchange implements the model of dual porosity.
IteratorBase end() const
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
EI_Address ei_address() const
Definition: accessors.cc:178
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:124
void setup_components()
MultiField< 3, FieldValue< 3 >::Scalar > isotherm_other
Langmuir sorption coeficients alpha (in fraction c_s = omega * (alpha*c_a)/(1- alpha*c_a)).
FieldFEScalarVec conc_solid_fe
Underlaying FieldFE for each substance of conc_solid.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:534
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
void output_data(void) override
Output method.
static constexpr Mask input_copy
Definition: field_flag.hh:44
MultiField< 3, FieldValue< 3 >::Scalar > init_conc_solid
Initial sorbed concentrations.
void initialize_from_input()
Initializes private members of sorption from the input record.
std::vector< double > solubility_vec_
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
unsigned int n_interpolation_steps_
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:93
void make_tables(void)
#define START_TIMER(tag)
Starts a timer with specified tag.
Mesh * mesh_
Definition: equation.hh:218
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.
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
void update_max_conc()
Computes maximal aqueous concentration at the current step.
virtual void compute_common_ele_data(const ElementAccessor< 3 > &elem)=0
MultiField< 3, FieldValue< 3 >::Enum > sorption_type
Discrete need Selection for initialization.
void update_solution(void) override
Updates the solution.
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:458
Region region() const
Definition: accessors.hh:165
SorptionType
Type of adsorption isotherm.
Definition: isotherm.hh:174
Definition: system.hh:65
Support classes for parallel programing.
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
Definition: field_set.hh:193
Class SorptionBase is abstract class representing model of sorption in transport. ...
void interpolate(double &c_aqua, double &c_sorbed)
Definition: isotherm.cc:84
void initialize(std::shared_ptr< OutputTime > stream, Mesh *mesh, Input::Record in_rec, const TimeGovernor &tg)
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
Definition: asserts.hh:336
const Selection & close() const
Close the Selection, no more values can be added.
Input::Record input_record_
Definition: equation.hh:220
void set_components(const std::vector< string > &names)
Definition: field_set.hh:180
void set_initial_condition()
Reads and sets initial condition for concentration in solid.
EqData(const string &output_field_name, const string &output_field_desc)
Collect all fields.
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than)
Definition: asserts.hh:296
FieldFEScalarVec conc_mobile_fe
FieldFEs representing P0 interpolation of mobile concentration (passed from transport).
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:157
std::shared_ptr< DOFHandlerMultiDim > dof_handler_
Pointer to DOF handler used through the reaction tree.
bool is_constant(Region reg) const
Definition: field_set.cc:173
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:91
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:258
#define END_TIMER(tag)
Ends a timer with specified tag.
void clear_max_conc()
Sets max conc to zeros on all regins.
bool is_precomputed(void)
Definition: isotherm.hh:239
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:186
void set_components(const std::vector< string > &names)
Record type proxy class.
Definition: type_record.hh:182
static const Input::Type::Record & get_input_type()
void isotherm_reinit_all(const ElementAccessor< 3 > &elm)
Calls isotherm_reinit for all isotherms.
Base of exceptions used in Flow123d.
Definition: exceptions.hh:75
MultiField< 3, FieldValue< 3 >::Scalar > distribution_coefficient
Multiplication coefficients (k, omega) for all types of isotherms.
unsigned int n_substances_
Class for representation SI units of Fields.
Definition: unit_si.hh:40
EqData * data_
Pointer to equation data. The object is constructed in descendants.
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:96
std::vector< std::vector< double > > max_conc
Other possible transformation of coordinates:
Class for declaration of the input data that are in string format.
Definition: type_base.hh:582
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:53
bool changed() const
Definition: field_set.cc:165
Template for classes storing finite set of named values.
Input::Type::Record make_field_descriptor_type(const std::string &equation_name) const
Definition: field_set.cc:62
struct SorptionBase::CommonElementData common_ele_data
TimeGovernor * time_
Definition: equation.hh:219
FieldSet input_data_set_
Input data set - fields in this set are read from the input file.
std::shared_ptr< ReactionTerm > reaction_liquid
void initialize() override
Prepares the object to usage.