Flow123d  release_2.2.0-914-gf1a3a4f
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 
26 #include "reaction/isotherm.hh"
27 
28 #include "system/system.hh"
29 #include "system/sys_profiler.hh"
30 
31 #include "la/distribution.hh"
32 #include "mesh/mesh.h"
33 #include "mesh/elements.h"
34 #include "input/type_selection.hh"
35 
36 #include "fields/field_set.hh"
38 
39 using namespace Input::Type;
40 
42  return Selection("SorptionType")
43  .add_value(Isotherm::none,"none", "No sorption considered.")
44  .add_value(Isotherm::linear, "linear",
45  "Linear isotherm runs the concentration exchange between liquid and solid.")
46  .add_value(Isotherm::langmuir, "langmuir",
47  "Langmuir isotherm runs the concentration exchange between liquid and solid.")
48  .add_value(Isotherm::freundlich, "freundlich",
49  "Freundlich isotherm runs the concentration exchange between liquid and solid.")
50  .close();
51 }
52 
53 
54 
56  return Record("SorptionAux", "AUXILIARY RECORD. Should not be directly part of the input tree.")
57  .declare_key("substances", Array(String(),1), Default::obligatory(),
58  "Names of the substances that take part in the sorption model.")
59  .declare_key("solvent_density", Double(0.0), Default("1.0"),
60  "Density of the solvent.")
61  .declare_key("substeps", Integer(1), Default("1000"),
62  "Number of equidistant substeps, molar mass and isotherm intersections")
63  .declare_key("solubility", Array(Double(0.0)), Default::optional(), //("-1.0"), //
64  "Specifies solubility limits of all the sorbing species.")
65  .declare_key("table_limits", Array(Double(0.0)), Default::optional(), //("-1.0"), //
66  "Specifies highest aqueous concentration in interpolation table.")
67  .declare_key("input_fields", Array(EqData("").input_data_set_.make_field_descriptor_type("Sorption")), Default::obligatory(), //
68  "Containes region specific data necessary to construct isotherms.")//;
69  .declare_key("reaction_liquid", ReactionTerm::it_abstract_reaction(), Default::optional(), "Reaction model following the sorption in the liquid.")
70  .declare_key("reaction_solid", ReactionTerm::it_abstract_reaction(), Default::optional(), "Reaction model following the sorption in the solid.")
71  .close();
72 }
73 
74 
75 SorptionBase::EqData::EqData(const string &output_field_name)
76 {
77  ADD_FIELD(rock_density, "Rock matrix density.", "0.0");
78  rock_density.units( UnitSI().kg().m(-3) );
79 
80  ADD_FIELD(sorption_type,"Considered sorption is described by selected isotherm. If porosity on an element is equal or even higher than 1.0 (meaning no sorbing surface), then type 'none' will be selected automatically."); //
81  sorption_type.input_selection(get_sorption_type_selection());
82  sorption_type.units( UnitSI::dimensionless() );
83 
84  ADD_FIELD(distribution_coefficient,"Multiplication parameters (k, omega) in either Langmuir c_s = omega * (alpha*c_a)/(1- alpha*c_a) or in linear c_s = k * c_a isothermal description.","1.0");
85  distribution_coefficient.units( UnitSI().m(3).kg(-1) );
86 
87  ADD_FIELD(isotherm_other,"Second parameters (alpha, ...) defining isotherm c_s = omega * (alpha*c_a)/(1- alpha*c_a).","1.0");
88  isotherm_other.units( UnitSI::dimensionless() );
89 
90  ADD_FIELD(init_conc_solid, "Initial solid concentration of substances."
91  " Vector, one value for every substance.", "0");
92  init_conc_solid.units( UnitSI().mol().kg(-1) );
93 
94  input_data_set_ += *this;
95 
96  // porosity field is set from governing equation (transport) later
97  // hence we do not add it to the input_data_set_
98  *this += porosity
99  .name("porosity")
100  .units( UnitSI::dimensionless() )
101  .flags(FieldFlag::input_copy)
102  .set_limits(0.0);
103 
104  output_fields += *this;
105  output_fields += conc_solid.name(output_field_name).units( UnitSI().dimensionless() ).flags(FieldFlag::equation_result);
106 }
107 
108 
110  : ReactionTerm(init_mesh, in_rec),
111  data_(nullptr)
112 {
113  // creating reaction from input and setting their parameters
114  make_reactions();
115 }
116 
117 
119 {
120  if (data_ != nullptr) delete data_;
121 
122  VecScatterDestroy(&(vconc_out_scatter));
123  if (vconc_solid != NULL) {
124 
125 
126  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
127  {
128  //no mpi vectors
129  VecDestroy( &(vconc_solid[sbi]) );
130  delete [] conc_solid[sbi];
131  }
132  delete [] vconc_solid;
133  delete [] conc_solid;
134  }
135 }
136 
138 {
140 
141  reactions_it = input_record_.find<Input::AbstractRecord>("reaction_liquid");
142  if ( reactions_it )
143  {
144  // TODO: allowed instances in this case are only
145  // FirstOrderReaction, RadioactiveDecay
146  reaction_liquid = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
147  } else
148  {
149  reaction_liquid = nullptr;
150  }
151 
152  reactions_it = input_record_.find<Input::AbstractRecord>("reaction_solid");
153  if ( reactions_it )
154  {
155  // TODO: allowed instances in this case are only
156  // FirstOrderReaction, RadioactiveDecay
157  reaction_solid = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
158  } else
159  {
160  reaction_solid = nullptr;
161  }
162 }
163 
165 {
166  OLD_ASSERT(distribution_ != nullptr, "Distribution has not been set yet.\n");
167  OLD_ASSERT(time_ != nullptr, "Time governor has not been set yet.\n");
168  OLD_ASSERT(output_stream_,"Null output stream.");
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  for(unsigned int i_reg = 0; i_reg < nr_of_regions; i_reg++)
178  {
179  isotherms[i_reg].resize(n_substances_);
180  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++)
181  {
182  isotherms[i_reg][i_subst] = Isotherm();
183  }
184  }
185 
186  //allocating new array for sorbed concentrations
187  conc_solid = new double* [substances_.size()];
188  conc_solid_out.clear();
189  conc_solid_out.resize( substances_.size() );
190  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
191  {
192  conc_solid[sbi] = new double [distribution_->lsize()];
193  conc_solid_out[sbi].resize( distribution_->size() );
194  //zero initialization of solid concentration for all substances
195  for(unsigned int i=0; i < distribution_->lsize(); i++)
196  conc_solid[sbi][i] = 0;
197  }
198 
200 
202 
203  if(reaction_liquid)
204  {
205  reaction_liquid->substances(substances_)
206  .concentration_matrix(concentration_matrix_, distribution_, el_4_loc_, row_4_el_)
207  .set_time_governor(*time_);
208  reaction_liquid->initialize();
209  }
210  if(reaction_solid)
211  {
212  reaction_solid->substances(substances_)
213  .concentration_matrix(conc_solid, distribution_, el_4_loc_, row_4_el_)
214  .set_time_governor(*time_);
215  reaction_solid->initialize();
216  }
217 }
218 
219 
221 {
222  Input::Array substances_array = input_record_.val<Input::Array>("substances");
223  unsigned int k, global_idx, i_subst = 0;
224  bool found;
225  Input::Iterator<string> spec_iter = substances_array.begin<string>();
226 
227  for(; spec_iter != substances_array.end(); ++spec_iter, i_subst++)
228  {
229  //finding the name of a substance in the global array of names
230  found = false;
231  for(k = 0; k < substances_.size(); k++)
232  {
233  if (*spec_iter == substances_[k].name())
234  {
235  global_idx = k;
236  found = true;
237  break;
238  }
239  }
240 
241  if(!found)
242  THROW(ReactionTerm::ExcUnknownSubstance()
243  << ReactionTerm::EI_Substance(*spec_iter)
244  << substances_array.ei_address());
245 
246  //finding the global index of substance in the local array
247  found = false;
248  for(k = 0; k < substance_global_idx_.size(); k++)
249  {
250  if(substance_global_idx_[k] == global_idx)
251  {
252  found = true;
253  break;
254  }
255  }
256 
257  if(!found)
258  {
259  substance_global_idx_.push_back(global_idx);
260  }
261 
262  }
264 }
265 
267 {
268  // read number of interpolation steps - value checked by the record definition
269  n_interpolation_steps_ = input_record_.val<int>("substeps");
270 
271  // read the density of solvent - value checked by the record definition
272  solvent_density_ = input_record_.val<double>("solvent_density");
273 
274  // read the solubility vector
275  Input::Iterator<Input::Array> solub_iter = input_record_.find<Input::Array>("solubility");
276  if( solub_iter )
277  {
278  if (solub_iter->Array::size() != n_substances_)
279  {
280  THROW(SorptionBase::ExcSubstanceCountMatch()
281  << SorptionBase::EI_ArrayName("solubility")
283  // there is no way to get ei_address from 'solub_iter', only as a string
284  }
285 
286  else solub_iter->copy_to(solubility_vec_);
287  }
288  else{
289  // fill solubility_vec_ with zeros
290  solubility_vec_.clear();
291  solubility_vec_.resize(n_substances_,0.0);
292  }
293 
294  // read the interpolation table limits
295  Input::Iterator<Input::Array> interp_table_limits = input_record_.find<Input::Array>("table_limits");
296  if( interp_table_limits )
297  {
298  if (interp_table_limits->Array::size() != n_substances_)
299  {
300  THROW(SorptionBase::ExcSubstanceCountMatch()
301  << SorptionBase::EI_ArrayName("table_limits")
303  // there is no way to get ei_address from 'interp_table_limits', only as a string
304  }
305 
306  else interp_table_limits->copy_to(table_limit_);
307  }
308  else{
309  // fill table_limit_ with zeros
310  table_limit_.clear();
311  table_limit_.resize(n_substances_,0.0);
312  }
313 }
314 
316 {
317  OLD_ASSERT(n_substances_ > 0, "Number of substances is wrong, they might have not been set yet.\n");
318 
319  // create vector of substances that are involved in sorption
320  // and initialize data_ with their names
321  std::vector<std::string> substances_sorption;
322  for (unsigned int i : substance_global_idx_)
323  substances_sorption.push_back(substances_[i].name());
324  data_->set_components(substances_sorption);
325 
326  // read fields from input file
328 
329  data_->set_mesh(*mesh_);
330 
331  //initialization of output
332  //output_array = input_record_.val<Input::Array>("output_fields");
337  for (unsigned int sbi=0; sbi<substances_.size(); sbi++)
338  {
339  // create shared pointer to a FieldElementwise and push this Field to output_field on all regions
340  auto output_field_ptr = conc_solid_out[sbi].create_field<3, FieldValue<3>::Scalar>(substances_.size());
341  data_->conc_solid[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
342  }
343  //output_stream_->add_admissible_field_names(output_array);
345 }
346 
347 
349 {
350  OLD_ASSERT(distribution_ != nullptr, "Distribution has not been set yet.\n");
351  OLD_ASSERT(time_ != nullptr, "Time governor has not been set yet.\n");
352  OLD_ASSERT(output_stream_,"Null output stream.");
354 
356  std::stringstream ss; // print warning message with table of uninitialized fields
357  if ( FieldCommon::print_message_table(ss, "sorption") ) {
358  WarningOut() << ss.str();
359  }
361  make_tables();
362 
363  // write initial condition
364  //output_vector_gather();
365  //data_->output_fields.set_time(time_->step(), LimitSide::right);
366  //data_->output_fields.output(output_stream_);
367 
368  if(reaction_liquid) reaction_liquid->zero_time_step();
369  if(reaction_solid) reaction_solid->zero_time_step();
370 
371  output_data();
372 }
373 
375 {
376  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
377  {
378  unsigned int index = el_4_loc_[loc_el];
379  ElementAccessor<3> ele_acc = mesh_->element_accessor(index);
380 
381  //setting initial solid concentration for substances involved in adsorption
382  for (unsigned int sbi = 0; sbi < n_substances_; sbi++)
383  {
384  int subst_id = substance_global_idx_[sbi];
385  conc_solid[subst_id][loc_el] = data_->init_conc_solid[sbi].value(ele_acc.centre(), ele_acc);
386  }
387  }
388 }
389 
390 
392 {
393  data_->set_time(time_->step(), LimitSide::right); // set to the last computed time
394 
395  // if parameters changed during last time step, reinit isotherms and eventualy
396  // update interpolation tables in the case of constant rock matrix parameters
397  if(data_->changed())
398  make_tables();
399 
400 
401  START_TIMER("Sorption");
402  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
403  {
405  }
406  END_TIMER("Sorption");
407 
408  if(reaction_liquid) reaction_liquid->update_solution();
409  if(reaction_solid) reaction_solid->update_solution();
410 }
411 
413 {
414  try
415  {
416  ElementAccessor<3> elm;
417  BOOST_FOREACH(const Region &reg_iter, this->mesh_->region_db().get_region_set("BULK") )
418  {
419  int reg_idx = reg_iter.bulk_idx();
420 
421  if(data_->is_constant(reg_iter))
422  {
423  ElementAccessor<3> elm(this->mesh_, reg_iter); // constant element accessor
424  isotherm_reinit(isotherms[reg_idx],elm);
425  for(unsigned int i_subst = 0; i_subst < n_substances_; i_subst++)
426  {
427  isotherms[reg_idx][i_subst].make_table(n_interpolation_steps_);
428  }
429  }
430  }
431  }
432  catch(ExceptionBase const &e)
433  {
434  e << input_record_.ei_address();
435  throw;
436  }
437 }
438 
439 double **SorptionBase::compute_reaction(double **concentrations, int loc_el)
440 {
441  ElementFullIter elem = mesh_->element(el_4_loc_[loc_el]);
442  int reg_idx = elem->region().bulk_idx();
443  unsigned int i_subst, subst_id;
444 
445  std::vector<Isotherm> & isotherms_vec = isotherms[reg_idx];
446 
447  try{
448  // Constant value of rock density and mobile porosity over the whole region
449  // => interpolation_table is precomputed
450  if (isotherms_vec[0].is_precomputed())
451  {
452  for(i_subst = 0; i_subst < n_substances_; i_subst++)
453  {
454  subst_id = substance_global_idx_[i_subst];
455 
456  isotherms_vec[i_subst].interpolate(concentration_matrix_[subst_id][loc_el],
457  conc_solid[subst_id][loc_el]);
458  }
459  }
460  else
461  {
462  isotherm_reinit(isotherms_vec, elem->element_accessor());
463 
464  for(i_subst = 0; i_subst < n_substances_; i_subst++)
465  {
466  subst_id = substance_global_idx_[i_subst];
467  isotherms_vec[i_subst].compute(concentration_matrix_[subst_id][loc_el],
468  conc_solid[subst_id][loc_el]);
469  }
470  }
471  }
472  catch(ExceptionBase const &e)
473  {
474  e << input_record_.ei_address();
475  throw;
476  }
477 
478  return concentrations;
479 }
480 
481 
482 /**************************************** OUTPUT ***************************************************/
483 
485 {
486  int sbi, n_subst;
487  n_subst = substances_.size();
488 
489  vconc_solid = new Vec [n_subst];
490 
491  for (sbi = 0; sbi < n_subst; sbi++) {
492  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, distribution_->lsize(), mesh_->n_elements(), conc_solid[sbi],
493  &vconc_solid[sbi]);
494  VecZeroEntries(vconc_solid[sbi]);
495 
496  VecZeroEntries(conc_solid_out[sbi].get_data_petsc());
497  }
498 
499  // creating output vector scatter
500  IS is;
501  ISCreateGeneral(PETSC_COMM_SELF, mesh_->n_elements(), row_4_el_, PETSC_COPY_VALUES, &is); //WithArray
502  VecScatterCreate(vconc_solid[0], is, conc_solid_out[0].get_data_petsc(), PETSC_NULL, &vconc_out_scatter);
503  ISDestroy(&(is));
504 }
505 
506 
508 {
509  unsigned int sbi;
510 
511  for (sbi = 0; sbi < substances_.size(); sbi++) {
512  VecScatterBegin(vconc_out_scatter, vconc_solid[sbi], conc_solid_out[sbi].get_data_petsc(), INSERT_VALUES, SCATTER_FORWARD);
513  VecScatterEnd(vconc_out_scatter, vconc_solid[sbi], conc_solid_out[sbi].get_data_petsc(), INSERT_VALUES, SCATTER_FORWARD);
514  }
515 }
516 
517 
519 {
523  }
524 
525  // Register fresh output data
526  data_->output_fields.output(time().step());
527 }
TimeGovernor & time()
Definition: equation.hh:148
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:201
unsigned int size() const
get global size
void set_input_list(Input::Array input_list)
Definition: field_set.hh:180
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.
void make_reactions()
std::vector< VectorSeqDouble > conc_solid_out
sorbed concentration array output (gathered - sequential)
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
double solvent_density_
EI_Address ei_address() const
Definition: accessors.cc:314
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: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
EqData(const string &output_field_name)
Collect all fields.
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
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
Definition: mesh.h:99
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:170
IdxInt * el_4_loc_
Indices of elements belonging to local dofs.
ElementAccessor< 3 > element_accessor(unsigned int idx, bool boundary=false)
Definition: mesh.cc:668
const TimeStep & step(int index=-1) const
std::vector< double > table_limit_
Class for declaration of the integral input data.
Definition: type_base.hh:489
static Input::Type::Abstract & it_abstract_reaction()
#define ADD_FIELD(name,...)
Definition: field_set.hh:269
unsigned int size() const
Definition: substance.hh:87
Class for declaration of inputs sequences.
Definition: type_base.hh:345
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.
Class Dual_por_exchange implements the model of dual porosity.
IteratorBase end() const
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
#define OLD_ASSERT(...)
Definition: global_defs.h:131
void setup_components()
unsigned int n_elements() const
Definition: mesh.h:156
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:540
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).
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
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
IdxInt * row_4_el_
Indices of rows belonging to elements.
Mesh * mesh_
Definition: equation.hh:223
virtual void isotherm_reinit(std::vector< Isotherm > &isotherms, const ElementAccessor< 3 > &elm)=0
Reinitializes the isotherm.
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:490
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
double ** compute_reaction(double **concentrations, int loc_el) override
Support classes for parallel programing.
Class SorptionBase is abstract class representing model of sorption in transport. ...
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:167
void set_initial_condition()
Reads and sets initial condition for concentration in solid.
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:90
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
#define END_TIMER(tag)
Ends a timer with specified tag.
arma::vec::fixed< spacedim > centre() const
Definition: accessors.hh:91
double ** conc_solid
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:173
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
Base of exceptions used in Flow123d.
Definition: exceptions.hh:75
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:91
Other possible transformation of coordinates:
Class for declaration of the input data that are in string format.
Definition: type_base.hh:588
#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
TimeGovernor * time_
Definition: equation.hh:224
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:260
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.