Flow123d  release_3.0.0-910-g56968e3
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 
18 #include <boost/foreach.hpp>
19 
25 #include "reaction/isotherm.hh"
26 
27 #include "system/system.hh"
28 #include "system/sys_profiler.hh"
29 
30 #include "la/distribution.hh"
31 #include "mesh/mesh.h"
32 #include "mesh/elements.h"
33 #include "mesh/accessors.hh"
34 #include "input/type_selection.hh"
35 
36 #include "fields/field_set.hh"
37 #include "fields/field_fe.hh"
39 
40 using namespace Input::Type;
41 
43  return Selection("SorptionType")
44  .add_value(Isotherm::none,"none", "No sorption considered.")
45  .add_value(Isotherm::linear, "linear",
46  "Linear isotherm runs the concentration exchange between liquid and solid.")
47  .add_value(Isotherm::langmuir, "langmuir",
48  "Langmuir isotherm runs the concentration exchange between liquid and solid.")
49  .add_value(Isotherm::freundlich, "freundlich",
50  "Freundlich isotherm runs the concentration exchange between liquid and solid.")
51  .close();
52 }
53 
54 
55 
57  return Record("SorptionAux", "AUXILIARY RECORD. Should not be directly part of the input tree.")
58  .declare_key("substances", Array(String(),1), Default::obligatory(),
59  "Names of the substances that take part in the sorption model.")
60  .declare_key("solvent_density", Double(0.0), Default("1.0"),
61  "Density of the solvent.")
62  .declare_key("substeps", Integer(1), Default("1000"),
63  "Number of equidistant substeps, molar mass and isotherm intersections")
64  .declare_key("solubility", Array(Double(0.0)), Default::optional(), //("-1.0"), //
65  "Specifies solubility limits of all the sorbing species.")
66  .declare_key("table_limits", Array(Double(-1.0)), Default::optional(), //("-1.0"), //
67  "Specifies the highest aqueous concentration in the isotherm function interpolation table. "
68  "Use any negative value for an automatic choice according to current maximal concentration (default and recommended). "
69  "Use '0' to always evaluate isotherm function directly (can be very slow). "
70  "Use a positive value to set the interpolation table limit manually "
71  "(if aqueous concentration is higher, then the isotherm function is evaluated directly).")
72  .declare_key("input_fields", Array(EqData("","").input_data_set_.make_field_descriptor_type("Sorption")), Default::obligatory(), //
73  "Containes region specific data necessary to construct isotherms.")//;
74  .declare_key("reaction_liquid", ReactionTerm::it_abstract_reaction(), Default::optional(), "Reaction model following the sorption in the liquid.")
75  .declare_key("reaction_solid", ReactionTerm::it_abstract_reaction(), Default::optional(), "Reaction model following the sorption in the solid.")
76  .close();
77 }
78 
79 
80 SorptionBase::EqData::EqData(const string &output_field_name, const string &output_field_desc)
81 {
82  *this += rock_density.name("rock_density")
83  .description("Rock matrix density.")
84  .input_default("0.0")
85  .units( UnitSI().kg().m(-3) );
86 
87  *this += sorption_type.name("sorption_type")
88  .description("Considered sorption is described by selected isotherm.\n"
89  "If porosity on an element is equal to 1.0 (or even higher), meaning no sorbing surface, then type 'none' will be selected automatically.")
90  .input_selection(get_sorption_type_selection())
91  .units( UnitSI::dimensionless() );
92 
93  *this += distribution_coefficient.name("distribution_coefficient")
94  .description("Distribution coefficient (( $k_l, k_F, k_L $)) of linear, Freundlich or Langmuir isotherm respectively.")
95  .input_default("1.0")
96  .units( UnitSI().m(3).kg(-1) );
97 
98  *this += isotherm_other.name("isotherm_other")
99  .description("Additional parameter (($ \\alpha $)) of nonlinear isotherms.")
100  .input_default("1.0")
101  .units( UnitSI::dimensionless() );
102 
103  *this += init_conc_solid.name("init_conc_solid")
104  .description("Initial solid concentration of substances. It is a vector: one value for every substance.")
105  .input_default("0")
106  .units( UnitSI().dimensionless() );
107 
108  input_data_set_ += *this;
109 
110  // porosity field is set from governing equation (transport) later
111  // hence we do not add it to the input_data_set_
112  *this += porosity
113  .name("porosity")
114  .units( UnitSI::dimensionless() )
115  .flags(FieldFlag::input_copy)
116  .set_limits(0.0);
117 
118  output_fields += *this;
119  output_fields += conc_solid.name(output_field_name)
120  .description(output_field_desc)
121  .units( UnitSI().dimensionless() )
123 }
124 
125 
127  : ReactionTerm(init_mesh, in_rec),
128  data_(nullptr)
129 {
130  // creating reaction from input and setting their parameters
131  make_reactions();
132 }
133 
134 
136 {
137  if (data_ != nullptr) delete data_;
138 
139  //for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
140  //{
141  // //no mpi vectors
142  // delete [] conc_solid[sbi];
143  //}
144  delete [] conc_solid;
145 }
146 
148 {
150 
151  reactions_it = input_record_.find<Input::AbstractRecord>("reaction_liquid");
152  if ( reactions_it )
153  {
154  // TODO: allowed instances in this case are only
155  // FirstOrderReaction, RadioactiveDecay
156  reaction_liquid = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
157  } else
158  {
159  reaction_liquid = nullptr;
160  }
161 
162  reactions_it = input_record_.find<Input::AbstractRecord>("reaction_solid");
163  if ( reactions_it )
164  {
165  // TODO: allowed instances in this case are only
166  // FirstOrderReaction, RadioactiveDecay
167  reaction_solid = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
168  } else
169  {
170  reaction_solid = nullptr;
171  }
172 }
173 
175 {
176  ASSERT_PTR(distribution_).error("Distribution has not been set yet.\n");
177  ASSERT_PTR(time_).error("Time governor has not been set yet.\n");
178  ASSERT(output_stream_).error("Null output stream.");
179  ASSERT_LT(0, substances_.size());
180 
181  initialize_substance_ids(); //computes present substances and sets indices
182  initialize_from_input(); //reads non-field data from input
183 
184  //isotherms array resized bellow
185  unsigned int nr_of_regions = mesh_->region_db().bulk_size();
186  isotherms.resize(nr_of_regions);
187  max_conc.resize(nr_of_regions);
188  for(unsigned int i_reg = 0; i_reg < nr_of_regions; i_reg++)
189  {
190  isotherms[i_reg].resize(n_substances_);
191  max_conc[i_reg].resize(n_substances_, 0.0);
192  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++)
193  {
194  isotherms[i_reg][i_subst] = Isotherm();
195  }
196  }
197 
198  //allocating new array for sorbed concentrations
199  conc_solid = new double* [substances_.size()];
200  conc_solid_out.clear();
201  conc_solid_out.resize( substances_.size() );
202  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
203  {
204  conc_solid[sbi] = new double [distribution_->lsize()];
205  //zero initialization of solid concentration for all substances
206  for(unsigned int i=0; i < distribution_->lsize(); i++)
207  conc_solid[sbi][i] = 0;
208  }
209 
211 
212  if(reaction_liquid)
213  {
214  reaction_liquid->substances(substances_)
215  .concentration_matrix(concentration_matrix_, distribution_, el_4_loc_, row_4_el_)
216  .set_dh(this->dof_handler_)
217  .set_time_governor(*time_);
218  reaction_liquid->initialize();
219  }
220  if(reaction_solid)
221  {
222  reaction_solid->substances(substances_)
223  .concentration_matrix(conc_solid, distribution_, el_4_loc_, row_4_el_)
224  .set_dh(this->dof_handler_)
225  .set_time_governor(*time_);
226  reaction_solid->initialize();
227  }
228 }
229 
230 
232 {
233  Input::Array substances_array = input_record_.val<Input::Array>("substances");
234  unsigned int k, global_idx, i_subst = 0;
235  bool found;
236  Input::Iterator<string> spec_iter = substances_array.begin<string>();
237 
238  for(; spec_iter != substances_array.end(); ++spec_iter, i_subst++)
239  {
240  //finding the name of a substance in the global array of names
241  found = false;
242  for(k = 0; k < substances_.size(); k++)
243  {
244  if (*spec_iter == substances_[k].name())
245  {
246  global_idx = k;
247  found = true;
248  break;
249  }
250  }
251 
252  if(!found)
253  THROW(ReactionTerm::ExcUnknownSubstance()
254  << ReactionTerm::EI_Substance(*spec_iter)
255  << substances_array.ei_address());
256 
257  //finding the global index of substance in the local array
258  found = false;
259  for(k = 0; k < substance_global_idx_.size(); k++)
260  {
261  if(substance_global_idx_[k] == global_idx)
262  {
263  found = true;
264  break;
265  }
266  }
267 
268  if(!found)
269  {
270  substance_global_idx_.push_back(global_idx);
271  }
272 
273  }
275 }
276 
278 {
279  // read number of interpolation steps - value checked by the record definition
280  n_interpolation_steps_ = input_record_.val<int>("substeps");
281 
282  // read the density of solvent - value checked by the record definition
283  solvent_density_ = input_record_.val<double>("solvent_density");
284 
285  // read the solubility vector
286  Input::Iterator<Input::Array> solub_iter = input_record_.find<Input::Array>("solubility");
287  if( solub_iter )
288  {
289  if (solub_iter->Array::size() != n_substances_)
290  {
291  THROW(SorptionBase::ExcSubstanceCountMatch()
292  << SorptionBase::EI_ArrayName("solubility")
294  // there is no way to get ei_address from 'solub_iter', only as a string
295  }
296 
297  else solub_iter->copy_to(solubility_vec_);
298  }
299  else{
300  // fill solubility_vec_ with zeros
301  solubility_vec_.clear();
302  solubility_vec_.resize(n_substances_,0.0);
303  }
304 
305  // read the interpolation table limits
306  Input::Iterator<Input::Array> interp_table_limits = input_record_.find<Input::Array>("table_limits");
307  if( interp_table_limits )
308  {
309  if (interp_table_limits->Array::size() != n_substances_)
310  {
311  THROW(SorptionBase::ExcSubstanceCountMatch()
312  << SorptionBase::EI_ArrayName("table_limits")
314  // there is no way to get ei_address from 'interp_table_limits', only as a string
315  }
316 
317  else interp_table_limits->copy_to(table_limit_);
318  }
319  else{
320  // fill table_limit_ with negative values -> automatic limit
321  table_limit_.clear();
322  table_limit_.resize(n_substances_,-1.0);
323  }
324 }
325 
327 {
328  ASSERT_GT(n_substances_, 0).error("Number of substances is wrong, they might have not been set yet.\n");
329 
330  // create vector of substances that are involved in sorption
331  // and initialize data_ with their names
332  std::vector<std::string> substances_sorption;
333  for (unsigned int i : substance_global_idx_)
334  substances_sorption.push_back(substances_[i].name());
335  data_->set_components(substances_sorption);
336 
337  // read fields from input file
339 
340  data_->set_mesh(*mesh_);
341 
342  //initialization of output
343  //output_array = input_record_.val<Input::Array>("output_fields");
348  for (unsigned int sbi=0; sbi<substances_.size(); sbi++)
349  {
350  // create shared pointer to a FieldFE and push this Field to output_field on all regions
351  std::shared_ptr<FieldFE<3, FieldValue<3>::Scalar> > output_field_ptr = make_shared< FieldFE<3, FieldValue<3>::Scalar> >();
352  conc_solid_out[sbi] = output_field_ptr->set_fe_data(this->dof_handler_);
353  data_->conc_solid[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
354  double *out_array;
355  VecGetArray(conc_solid_out[sbi].petsc_vec(), &out_array);
356  conc_solid[sbi] = out_array;
357  VecRestoreArray(conc_solid_out[sbi].petsc_vec(), &out_array);
358  }
359  //output_stream_->add_admissible_field_names(output_array);
361 }
362 
363 
365 {
366  ASSERT_PTR(distribution_).error("Distribution has not been set yet.\n");
367  ASSERT_PTR(time_).error("Time governor has not been set yet.\n");
368  ASSERT(output_stream_).error("Null output stream.");
369  ASSERT_LT(0, substances_.size());
370 
372  std::stringstream ss; // print warning message with table of uninitialized fields
373  if ( FieldCommon::print_message_table(ss, "sorption") ) {
374  WarningOut() << ss.str();
375  }
377 
378  update_max_conc();
379  make_tables();
380 
381  // write initial condition
382  //data_->output_fields.set_time(time_->step(), LimitSide::right);
383  //data_->output_fields.output(output_stream_);
384 
385  if(reaction_liquid) reaction_liquid->zero_time_step();
386  if(reaction_solid) reaction_solid->zero_time_step();
387 
388  output_data();
389 }
390 
392 {
393  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
394  {
395  unsigned int index = el_4_loc_[loc_el];
396  ElementAccessor<3> ele_acc = mesh_->element_accessor(index);
397 
398  //setting initial solid concentration for substances involved in adsorption
399  for (unsigned int sbi = 0; sbi < n_substances_; sbi++)
400  {
401  int subst_id = substance_global_idx_[sbi];
402  conc_solid[subst_id][loc_el] = data_->init_conc_solid[sbi].value(ele_acc.centre(), ele_acc);
403  }
404  }
405 }
406 
407 
409 {
410  data_->set_time(time_->step(), LimitSide::right); // set to the last computed time
411 
412  // if parameters changed during last time step, reinit isotherms and eventualy
413  // update interpolation tables in the case of constant rock matrix parameters
414  make_tables();
415  clear_max_conc();
416 
417  START_TIMER("Sorption");
418  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
419  {
421  }
422  END_TIMER("Sorption");
423 
424  if(reaction_liquid) reaction_liquid->update_solution();
425  if(reaction_solid) reaction_solid->update_solution();
426 }
427 
428 void SorptionBase::isotherm_reinit(unsigned int i_subst, const ElementAccessor<3> &elem)
429 {
430  START_TIMER("SorptionBase::isotherm_reinit");
431 
432  double mult_coef = data_->distribution_coefficient[i_subst].value(elem.centre(),elem);
433  double second_coef = data_->isotherm_other[i_subst].value(elem.centre(),elem);
434 
435  int reg_idx = elem.region().bulk_idx();
436  Isotherm & isotherm = isotherms[reg_idx][i_subst];
437 
438  bool limited_solubility_on = solubility_vec_[i_subst] > 0.0;
439 
440  // in case of no sorbing surface, set isotherm type None
441  if( common_ele_data.no_sorbing_surface_cond <= std::numeric_limits<double>::epsilon())
442  {
443  isotherm.reinit(Isotherm::none, false, solvent_density_,
445  0,0,0);
446  return;
447  }
448 
449  if ( common_ele_data.scale_sorbed <= 0.0)
450  xprintf(UsrErr, "Scaling parameter in sorption is not positive. Check the input for rock density and molar mass of %d. substance.",i_subst);
451 
452  isotherm.reinit(Isotherm::SorptionType(data_->sorption_type[i_subst].value(elem.centre(),elem)),
453  limited_solubility_on, solvent_density_,
455  solubility_vec_[i_subst], mult_coef, second_coef);
456 }
457 
459 {
460  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++)
461  {
462  isotherm_reinit(i_subst, elem);
463  }
464 }
465 
467 {
468  unsigned int reg_idx, i_subst;
469 
470  // clear max concetrations array
471  unsigned int nr_of_regions = mesh_->region_db().bulk_size();
472  for(reg_idx = 0; reg_idx < nr_of_regions; reg_idx++)
473  for(i_subst = 0; i_subst < n_substances_; i_subst++)
474  max_conc[reg_idx][i_subst] = 0.0;
475 }
476 
478 {
479  unsigned int reg_idx, i_subst, subst_id;
480 
481  clear_max_conc();
482 
483  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++){
485  reg_idx = ele.region().bulk_idx();
486  for(i_subst = 0; i_subst < n_substances_; i_subst++){
487  subst_id = substance_global_idx_[i_subst];
488  max_conc[reg_idx][i_subst] = std::max(max_conc[reg_idx][i_subst], concentration_matrix_[subst_id][loc_el]);
489  }
490  }
491 }
492 
494 {
495  START_TIMER("SorptionBase::make_tables");
496  try
497  {
498  ElementAccessor<3> elm;
499  for(const Region &reg_iter: this->mesh_->region_db().get_region_set("BULK"))
500  {
501  int reg_idx = reg_iter.bulk_idx();
502  // true if data has been changed and are constant on the region
503  bool call_reinit = data_->changed() && data_->is_constant(reg_iter);
504 
505  if(call_reinit)
506  {
507  ElementAccessor<3> elm(this->mesh_, reg_iter);
508 // DebugOut().fmt("isotherm reinit\n");
510  isotherm_reinit_all(elm);
511  }
512 
513  // find table limit and create interpolation table for every substance
514  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++){
515 
516  // clear interpolation tables, if not spacially constant OR switched off
517  if(! data_->is_constant(reg_iter) || table_limit_[i_subst] == 0.0){
518  isotherms[reg_idx][i_subst].clear_table();
519 // DebugOut().fmt("limit: 0.0 -> clear table\n");
520  continue;
521  }
522 
523  // if true then make_table will be called at the end
524  bool call_make_table = call_reinit;
525  // initialy try to keep the current table limit (it is zero at zero time step)
526  double subst_table_limit = isotherms[reg_idx][i_subst].table_limit();
527 
528  // if automatic, possibly remake tables with doubled range when table maximum was reached
529  if(table_limit_[i_subst] < 0.0)
530  {
531  if(subst_table_limit < max_conc[reg_idx][i_subst])
532  {
533  call_make_table = true;
534  subst_table_limit = 2*max_conc[reg_idx][i_subst];
535 // DebugOut().fmt("limit: max conc\n");
536  }
537  }
538  // if not automatic, set given table limit
539  else
540  {
541  subst_table_limit = table_limit_[i_subst];
542  }
543 
544  if(call_make_table){
545  isotherms[reg_idx][i_subst].make_table(n_interpolation_steps_, subst_table_limit);
546 // DebugOut().fmt("reg: {} i_subst {}: table_limit = {}\n", reg_idx, i_subst, isotherms[reg_idx][i_subst].table_limit());
547  }
548  }
549  }
550  }
551  catch(ExceptionBase const &e)
552  {
553  e << input_record_.ei_address();
554  throw;
555  }
556 }
557 
558 double **SorptionBase::compute_reaction(double **concentrations, int loc_el)
559 {
561  int reg_idx = elem.region().bulk_idx();
562  unsigned int i_subst, subst_id;
563  // for checking, that the common element data are computed once at maximum
564  bool is_common_ele_data_valid = false;
565 
566  try{
567  for(i_subst = 0; i_subst < n_substances_; i_subst++)
568  {
569  subst_id = substance_global_idx_[i_subst];
570  Isotherm & isotherm = isotherms[reg_idx][i_subst];
571  if (isotherm.is_precomputed()){
572 // DebugOut().fmt("isotherms precomputed - interpolate, subst[{}]\n", i_subst);
573  isotherm.interpolate(concentration_matrix_[subst_id][loc_el],
574  conc_solid[subst_id][loc_el]);
575  }
576  else{
577 // DebugOut().fmt("isotherms reinit - compute , subst[{}]\n", i_subst);
578  if(! is_common_ele_data_valid){
580  is_common_ele_data_valid = true;
581  }
582 
583  isotherm_reinit(i_subst, elem);
584  isotherm.compute(concentration_matrix_[subst_id][loc_el],
585  conc_solid[subst_id][loc_el]);
586  }
587 
588  // update maximal concentration per region (optimization for interpolation)
589  if(table_limit_[i_subst] < 0)
590  max_conc[reg_idx][i_subst] = std::max(max_conc[reg_idx][i_subst],
591  concentration_matrix_[subst_id][loc_el]);
592  }
593  }
594  catch(ExceptionBase const &e)
595  {
596  e << input_record_.ei_address();
597  throw;
598  }
599 
600  return concentrations;
601 }
602 
603 
604 /**************************************** OUTPUT ***************************************************/
605 
607 {
609  // Register fresh output data
610  data_->output_fields.output(time().step());
611 }
TimeGovernor & time()
Definition: equation.hh:148
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:211
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:329
void make_reactions()
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
double solvent_density_
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.
double ** concentration_matrix_
unsigned int bulk_size() const
Definition: region.cc:269
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
void output(TimeStep step)
#define ASSERT_GT(a, b)
Definition of comparative assert macro (Greater Than)
Definition: asserts.hh:311
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:76
Iterator< Ret > find(const string &key) const
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()
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:346
const TimeStep & step(int index=-1) const
std::vector< double > table_limit_
Class for declaration of the integral input data.
Definition: type_base.hh:490
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:346
virtual ElementAccessor< 3 > element_accessor(unsigned int idx) const
Create and return ElementAccessor to element of given idx.
Definition: mesh.cc:726
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.
Definition: accessors.hh:285
EI_Address ei_address() const
Definition: accessors.cc:178
std::vector< VectorMPI > conc_solid_out
sorbed concentration array output (gathered - sequential)
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)).
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:541
Distribution * distribution_
Pointer to reference to distribution of elements between processors.
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
std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > output_field_ptr
Fields correspond with conc_solid_out.
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:292
unsigned int n_interpolation_steps_
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:92
void make_tables(void)
#define START_TIMER(tag)
Starts a timer with specified tag.
Mesh * mesh_
Definition: equation.hh:223
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:501
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:459
Region region() const
Definition: accessors.hh:95
SorptionType
Type of adsorption isotherm.
Definition: isotherm.hh:174
double ** compute_reaction(double **concentrations, int loc_el) override
Definition: system.hh:64
Support classes for parallel programing.
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
Definition: field_set.hh:190
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:335
const Selection & close() const
Close the Selection, no more values can be added.
Input::Record input_record_
Definition: equation.hh:225
void set_components(const std::vector< string > &names)
Definition: field_set.hh:177
void set_initial_condition()
Reads and sets initial condition for concentration in solid.
LongIdx * row_4_el_
Indices of rows belonging to elements.
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:295
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:246
#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
double ** conc_solid
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:183
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.
LongIdx * el_4_loc_
Indices of elements belonging to local dofs.
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:589
#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:224
FieldSet input_data_set_
Input data set - fields in this set are read from the input file.
std::shared_ptr< ReactionTerm > reaction_liquid
unsigned int lsize(int proc) const
get local size
void initialize() override
Prepares the object to usage.