Flow123d  release_3.0.0-680-gbed5aba
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  chkerr(VecScatterDestroy(&(vconc_out_scatter)));
140  if (vconc_solid != NULL) {
141 
142 
143  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
144  {
145  //no mpi vectors
146  chkerr(VecDestroy( &(vconc_solid[sbi]) ));
147  delete [] conc_solid[sbi];
148  }
149  delete [] vconc_solid;
150  delete [] conc_solid;
151  }
152 }
153 
155 {
157 
158  reactions_it = input_record_.find<Input::AbstractRecord>("reaction_liquid");
159  if ( reactions_it )
160  {
161  // TODO: allowed instances in this case are only
162  // FirstOrderReaction, RadioactiveDecay
163  reaction_liquid = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
164  } else
165  {
166  reaction_liquid = nullptr;
167  }
168 
169  reactions_it = input_record_.find<Input::AbstractRecord>("reaction_solid");
170  if ( reactions_it )
171  {
172  // TODO: allowed instances in this case are only
173  // FirstOrderReaction, RadioactiveDecay
174  reaction_solid = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
175  } else
176  {
177  reaction_solid = nullptr;
178  }
179 }
180 
182 {
183  OLD_ASSERT(distribution_ != nullptr, "Distribution has not been set yet.\n");
184  OLD_ASSERT(time_ != nullptr, "Time governor has not been set yet.\n");
185  OLD_ASSERT(output_stream_,"Null output stream.");
187 
188  initialize_substance_ids(); //computes present substances and sets indices
189  initialize_from_input(); //reads non-field data from input
190 
191  //isotherms array resized bellow
192  unsigned int nr_of_regions = mesh_->region_db().bulk_size();
193  isotherms.resize(nr_of_regions);
194  max_conc.resize(nr_of_regions);
195  for(unsigned int i_reg = 0; i_reg < nr_of_regions; i_reg++)
196  {
197  isotherms[i_reg].resize(n_substances_);
198  max_conc[i_reg].resize(n_substances_, 0.0);
199  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++)
200  {
201  isotherms[i_reg][i_subst] = Isotherm();
202  }
203  }
204 
205  //allocating new array for sorbed concentrations
206  conc_solid = new double* [substances_.size()];
207  conc_solid_out.clear();
208  conc_solid_out.resize( substances_.size() );
209  output_field_ptr.clear();
210  output_field_ptr.resize( substances_.size() );
211  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
212  {
213  conc_solid[sbi] = new double [distribution_->lsize()];
214  conc_solid_out[sbi].resize( distribution_->size() );
215  //zero initialization of solid concentration for all substances
216  for(unsigned int i=0; i < distribution_->lsize(); i++)
217  conc_solid[sbi][i] = 0;
218  }
219 
221 
223 
224  if(reaction_liquid)
225  {
226  reaction_liquid->substances(substances_)
227  .concentration_matrix(concentration_matrix_, distribution_, el_4_loc_, row_4_el_)
228  .set_time_governor(*time_);
229  reaction_liquid->initialize();
230  }
231  if(reaction_solid)
232  {
233  reaction_solid->substances(substances_)
234  .concentration_matrix(conc_solid, distribution_, el_4_loc_, row_4_el_)
235  .set_time_governor(*time_);
236  reaction_solid->initialize();
237  }
238 }
239 
240 
242 {
243  Input::Array substances_array = input_record_.val<Input::Array>("substances");
244  unsigned int k, global_idx, i_subst = 0;
245  bool found;
246  Input::Iterator<string> spec_iter = substances_array.begin<string>();
247 
248  for(; spec_iter != substances_array.end(); ++spec_iter, i_subst++)
249  {
250  //finding the name of a substance in the global array of names
251  found = false;
252  for(k = 0; k < substances_.size(); k++)
253  {
254  if (*spec_iter == substances_[k].name())
255  {
256  global_idx = k;
257  found = true;
258  break;
259  }
260  }
261 
262  if(!found)
263  THROW(ReactionTerm::ExcUnknownSubstance()
264  << ReactionTerm::EI_Substance(*spec_iter)
265  << substances_array.ei_address());
266 
267  //finding the global index of substance in the local array
268  found = false;
269  for(k = 0; k < substance_global_idx_.size(); k++)
270  {
271  if(substance_global_idx_[k] == global_idx)
272  {
273  found = true;
274  break;
275  }
276  }
277 
278  if(!found)
279  {
280  substance_global_idx_.push_back(global_idx);
281  }
282 
283  }
285 }
286 
288 {
289  // read number of interpolation steps - value checked by the record definition
290  n_interpolation_steps_ = input_record_.val<int>("substeps");
291 
292  // read the density of solvent - value checked by the record definition
293  solvent_density_ = input_record_.val<double>("solvent_density");
294 
295  // read the solubility vector
296  Input::Iterator<Input::Array> solub_iter = input_record_.find<Input::Array>("solubility");
297  if( solub_iter )
298  {
299  if (solub_iter->Array::size() != n_substances_)
300  {
301  THROW(SorptionBase::ExcSubstanceCountMatch()
302  << SorptionBase::EI_ArrayName("solubility")
304  // there is no way to get ei_address from 'solub_iter', only as a string
305  }
306 
307  else solub_iter->copy_to(solubility_vec_);
308  }
309  else{
310  // fill solubility_vec_ with zeros
311  solubility_vec_.clear();
312  solubility_vec_.resize(n_substances_,0.0);
313  }
314 
315  // read the interpolation table limits
316  Input::Iterator<Input::Array> interp_table_limits = input_record_.find<Input::Array>("table_limits");
317  if( interp_table_limits )
318  {
319  if (interp_table_limits->Array::size() != n_substances_)
320  {
321  THROW(SorptionBase::ExcSubstanceCountMatch()
322  << SorptionBase::EI_ArrayName("table_limits")
324  // there is no way to get ei_address from 'interp_table_limits', only as a string
325  }
326 
327  else interp_table_limits->copy_to(table_limit_);
328  }
329  else{
330  // fill table_limit_ with negative values -> automatic limit
331  table_limit_.clear();
332  table_limit_.resize(n_substances_,-1.0);
333  }
334 }
335 
337 {
338  OLD_ASSERT(n_substances_ > 0, "Number of substances is wrong, they might have not been set yet.\n");
339 
340  // create vector of substances that are involved in sorption
341  // and initialize data_ with their names
342  std::vector<std::string> substances_sorption;
343  for (unsigned int i : substance_global_idx_)
344  substances_sorption.push_back(substances_[i].name());
345  data_->set_components(substances_sorption);
346 
347  // read fields from input file
349 
350  data_->set_mesh(*mesh_);
351 
352  //initialization of output
353  //output_array = input_record_.val<Input::Array>("output_fields");
358  for (unsigned int sbi=0; sbi<substances_.size(); sbi++)
359  {
360  // create shared pointer to a FieldFE and push this Field to output_field on all regions
361  output_field_ptr[sbi] = create_field<3, FieldValue<3>::Scalar>(conc_solid_out[sbi], *mesh_, 1);
362  data_->conc_solid[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr[sbi], 0);
363  }
364  //output_stream_->add_admissible_field_names(output_array);
366 }
367 
368 
370 {
371  OLD_ASSERT(distribution_ != nullptr, "Distribution has not been set yet.\n");
372  OLD_ASSERT(time_ != nullptr, "Time governor has not been set yet.\n");
373  OLD_ASSERT(output_stream_,"Null output stream.");
375 
377  std::stringstream ss; // print warning message with table of uninitialized fields
378  if ( FieldCommon::print_message_table(ss, "sorption") ) {
379  WarningOut() << ss.str();
380  }
382 
383  update_max_conc();
384  make_tables();
385 
386  // write initial condition
387  //output_vector_gather();
388  //data_->output_fields.set_time(time_->step(), LimitSide::right);
389  //data_->output_fields.output(output_stream_);
390 
391  if(reaction_liquid) reaction_liquid->zero_time_step();
392  if(reaction_solid) reaction_solid->zero_time_step();
393 
394  output_data();
395 }
396 
398 {
399  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
400  {
401  unsigned int index = el_4_loc_[loc_el];
402  ElementAccessor<3> ele_acc = mesh_->element_accessor(index);
403 
404  //setting initial solid concentration for substances involved in adsorption
405  for (unsigned int sbi = 0; sbi < n_substances_; sbi++)
406  {
407  int subst_id = substance_global_idx_[sbi];
408  conc_solid[subst_id][loc_el] = data_->init_conc_solid[sbi].value(ele_acc.centre(), ele_acc);
409  }
410  }
411 }
412 
413 
415 {
416  data_->set_time(time_->step(), LimitSide::right); // set to the last computed time
417 
418  // if parameters changed during last time step, reinit isotherms and eventualy
419  // update interpolation tables in the case of constant rock matrix parameters
420  make_tables();
421  clear_max_conc();
422 
423  START_TIMER("Sorption");
424  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
425  {
427  }
428  END_TIMER("Sorption");
429 
430  if(reaction_liquid) reaction_liquid->update_solution();
431  if(reaction_solid) reaction_solid->update_solution();
432 }
433 
434 void SorptionBase::isotherm_reinit(unsigned int i_subst, const ElementAccessor<3> &elem)
435 {
436  START_TIMER("SorptionBase::isotherm_reinit");
437 
438  double mult_coef = data_->distribution_coefficient[i_subst].value(elem.centre(),elem);
439  double second_coef = data_->isotherm_other[i_subst].value(elem.centre(),elem);
440 
441  int reg_idx = elem.region().bulk_idx();
442  Isotherm & isotherm = isotherms[reg_idx][i_subst];
443 
444  bool limited_solubility_on = solubility_vec_[i_subst] > 0.0;
445 
446  // in case of no sorbing surface, set isotherm type None
447  if( common_ele_data.no_sorbing_surface_cond <= std::numeric_limits<double>::epsilon())
448  {
449  isotherm.reinit(Isotherm::none, false, solvent_density_,
451  0,0,0);
452  return;
453  }
454 
455  if ( common_ele_data.scale_sorbed <= 0.0)
456  xprintf(UsrErr, "Scaling parameter in sorption is not positive. Check the input for rock density and molar mass of %d. substance.",i_subst);
457 
458  isotherm.reinit(Isotherm::SorptionType(data_->sorption_type[i_subst].value(elem.centre(),elem)),
459  limited_solubility_on, solvent_density_,
461  solubility_vec_[i_subst], mult_coef, second_coef);
462 }
463 
465 {
466  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++)
467  {
468  isotherm_reinit(i_subst, elem);
469  }
470 }
471 
473 {
474  unsigned int reg_idx, i_subst, subst_id;
475 
476  // clear max concetrations array
477  unsigned int nr_of_regions = mesh_->region_db().bulk_size();
478  for(reg_idx = 0; reg_idx < nr_of_regions; reg_idx++)
479  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++)
480  max_conc[reg_idx][i_subst] = 0.0;
481 }
482 
484 {
485  unsigned int reg_idx, i_subst, subst_id;
486 
487  clear_max_conc();
488 
489  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++){
491  reg_idx = ele.region().bulk_idx();
492  for(i_subst = 0; i_subst < n_substances_; i_subst++){
493  subst_id = substance_global_idx_[i_subst];
494  max_conc[reg_idx][i_subst] = std::max(max_conc[reg_idx][i_subst], concentration_matrix_[subst_id][loc_el]);
495  }
496  }
497 }
498 
500 {
501  START_TIMER("SorptionBase::make_tables");
502  try
503  {
504  ElementAccessor<3> elm;
505  for(const Region &reg_iter: this->mesh_->region_db().get_region_set("BULK"))
506  {
507  int reg_idx = reg_iter.bulk_idx();
508  // true if data has been changed and are constant on the region
509  bool call_reinit = data_->changed() && data_->is_constant(reg_iter);
510 
511  if(call_reinit)
512  {
513  ElementAccessor<3> elm(this->mesh_, reg_iter);
514 // DebugOut().fmt("isotherm reinit\n");
516  isotherm_reinit_all(elm);
517  }
518 
519  // find table limit and create interpolation table for every substance
520  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++){
521 
522  // clear interpolation tables, if not spacially constant OR switched off
523  if(! data_->is_constant(reg_iter) || table_limit_[i_subst] == 0.0){
524  isotherms[reg_idx][i_subst].clear_table();
525 // DebugOut().fmt("limit: 0.0 -> clear table\n");
526  continue;
527  }
528 
529  // if true then make_table will be called at the end
530  bool call_make_table = call_reinit;
531  // initialy try to keep the current table limit (it is zero at zero time step)
532  double subst_table_limit = isotherms[reg_idx][i_subst].table_limit();
533 
534  // if automatic, possibly remake tables with doubled range when table maximum was reached
535  if(table_limit_[i_subst] < 0.0)
536  {
537  if(subst_table_limit < max_conc[reg_idx][i_subst])
538  {
539  call_make_table = true;
540  subst_table_limit = 2*max_conc[reg_idx][i_subst];
541 // DebugOut().fmt("limit: max conc\n");
542  }
543  }
544  // if not automatic, set given table limit
545  else
546  {
547  subst_table_limit = table_limit_[i_subst];
548  }
549 
550  if(call_make_table){
551  isotherms[reg_idx][i_subst].make_table(n_interpolation_steps_, subst_table_limit);
552 // DebugOut().fmt("reg: {} i_subst {}: table_limit = {}\n", reg_idx, i_subst, isotherms[reg_idx][i_subst].table_limit());
553  }
554  }
555  }
556  }
557  catch(ExceptionBase const &e)
558  {
559  e << input_record_.ei_address();
560  throw;
561  }
562 }
563 
564 double **SorptionBase::compute_reaction(double **concentrations, int loc_el)
565 {
567  int reg_idx = elem.region().bulk_idx();
568  unsigned int i_subst, subst_id;
569  // for checking, that the common element data are computed once at maximum
570  bool is_common_ele_data_valid = false;
571 
572  try{
573  for(i_subst = 0; i_subst < n_substances_; i_subst++)
574  {
575  subst_id = substance_global_idx_[i_subst];
576  Isotherm & isotherm = isotherms[reg_idx][i_subst];
577  if (isotherm.is_precomputed()){
578 // DebugOut().fmt("isotherms precomputed - interpolate, subst[{}]\n", i_subst);
579  isotherm.interpolate(concentration_matrix_[subst_id][loc_el],
580  conc_solid[subst_id][loc_el]);
581  }
582  else{
583 // DebugOut().fmt("isotherms reinit - compute , subst[{}]\n", i_subst);
584  if(! is_common_ele_data_valid){
586  is_common_ele_data_valid = true;
587  }
588 
589  isotherm_reinit(i_subst, elem);
590  isotherm.compute(concentration_matrix_[subst_id][loc_el],
591  conc_solid[subst_id][loc_el]);
592  }
593 
594  // update maximal concentration per region (optimization for interpolation)
595  if(table_limit_[i_subst] < 0)
596  max_conc[reg_idx][i_subst] = std::max(max_conc[reg_idx][i_subst],
597  concentration_matrix_[subst_id][loc_el]);
598  }
599  }
600  catch(ExceptionBase const &e)
601  {
602  e << input_record_.ei_address();
603  throw;
604  }
605 
606  return concentrations;
607 }
608 
609 
610 /**************************************** OUTPUT ***************************************************/
611 
613 {
614  int sbi, n_subst;
615  n_subst = substances_.size();
616 
617  vconc_solid = new Vec [n_subst];
618 
619  for (sbi = 0; sbi < n_subst; sbi++) {
620  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, distribution_->lsize(), mesh_->n_elements(), conc_solid[sbi],
621  &vconc_solid[sbi]);
622  VecZeroEntries(vconc_solid[sbi]);
623 
624  VecZeroEntries(conc_solid_out[sbi].petsc_vec());
625  }
626 
627  // creating output vector scatter
628  IS is;
629  ISCreateGeneral(PETSC_COMM_SELF, mesh_->n_elements(), row_4_el_, PETSC_COPY_VALUES, &is); //WithArray
630  VecScatterCreate(vconc_solid[0], is, conc_solid_out[0].petsc_vec(), PETSC_NULL, &vconc_out_scatter);
631  ISDestroy(&(is));
632 }
633 
634 
636 {
637  unsigned int sbi;
638 
639  for (sbi = 0; sbi < substances_.size(); sbi++) {
640  VecScatterBegin(vconc_out_scatter, vconc_solid[sbi], conc_solid_out[sbi].petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
641  VecScatterEnd(vconc_out_scatter, vconc_solid[sbi], conc_solid_out[sbi].petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
642  }
643 }
644 
645 
647 {
651  }
652 
653  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++) {
655  }
656 
657  // Register fresh output data
658  data_->output_fields.output(time().step());
659 }
TimeGovernor & time()
Definition: equation.hh:148
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:211
unsigned int size() const
get global size
std::shared_ptr< ReactionTerm > reaction_solid
Iterator< ValueType > begin() const
void allocate_output_mpi(void)
Allocates petsc vectors, prepares them for output and creates vector scatter.
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)
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:80
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()
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:147
void fill_output_data(VectorMPI &vec_seq, std::shared_ptr< FieldFE< spacedim, Value > > field_ptr)
Definition: field_fe.hh:330
const RegionDB & region_db() const
Definition: mesh.h:147
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:719
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
#define OLD_ASSERT(...)
Definition: global_defs.h:131
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
Vec * vconc_solid
PETSC sorbed concentration vector (parallel).
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.
#define OLD_ASSERT_LESS(a, b)
Definition: global_defs.h:134
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 output_vector_gather(void) override
Gathers all the parallel vectors to enable them to be output.
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
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
Definition: mesh.h:346
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)
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.
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:157
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()
VecScatter vconc_out_scatter
Output vector scatter.
bool is_field_output_time(const FieldCommon &field, TimeStep step) const
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.