Flow123d
sorption_base.cc
Go to the documentation of this file.
1 #include <boost/foreach.hpp>
2 
3 #include "reaction/reaction.hh"
8 #include "reaction/isotherm.hh"
9 #include "reaction/sorption.hh"
10 
11 #include "system/system.hh"
12 #include "system/sys_profiler.hh"
13 
14 #include "la/distribution.hh"
15 #include "mesh/mesh.h"
16 #include "mesh/elements.h"
17 #include "input/type_selection.hh"
18 
19 #include "fields/field_set.hh"
21 
22 
23 using namespace std;
24 using namespace Input::Type;
25 
27  .add_value(Isotherm::none,"none", "No adsorption considered.")
28  .add_value(Isotherm::linear, "linear",
29  "Linear isotherm runs the concentration exchange between liquid and solid.")
30  .add_value(Isotherm::langmuir, "langmuir",
31  "Langmuir isotherm runs the concentration exchange between liquid and solid.")
32  .add_value(Isotherm::freundlich, "freundlich",
33  "Freundlich isotherm runs the concentration exchange between liquid and solid.");
34 
35 
36 
38  = Record("Adsorption", "AUXILIARY RECORD. Should not be directly part of the input tree.")
39  .declare_key("substances", Array(String()), Default::obligatory(),
40  "Names of the substances that take part in the adsorption model.")
41  .declare_key("solvent_density", Double(), Default("1.0"),
42  "Density of the solvent.")
43  .declare_key("substeps", Integer(), Default("1000"),
44  "Number of equidistant substeps, molar mass and isotherm intersections")
45  .declare_key("molar_mass", Array(Double()), Default::obligatory(),
46  "Specifies molar masses of all the adsorbing species.")
47  .declare_key("solubility", Array(Double(0.0)), Default::optional(), //("-1.0"), //
48  "Specifies solubility limits of all the adsorbing species.")
49  .declare_key("table_limits", Array(Double(0.0)), Default::optional(), //("-1.0"), //
50  "Specifies highest aqueous concentration in interpolation table.")
51  .declare_key("input_fields", Array(EqData("").input_data_set_.make_field_descriptor_type("Sorption")), Default::obligatory(), //
52  "Containes region specific data necessary to construct isotherms.")//;
53  .declare_key("reaction", ReactionTerm::input_type, Default::optional(), "Reaction model following the sorption.");
54 
55 
56 SorptionBase::EqData::EqData(const string &output_field_name)
57 {
58  ADD_FIELD(rock_density, "Rock matrix density.", "0.0");
59 
60  ADD_FIELD(sorption_type,"Considered adsorption is described by selected isotherm."); //
61  sorption_type.input_selection(&sorption_type_selection);
62 
63  ADD_FIELD(isotherm_mult,"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");
64 
65  ADD_FIELD(isotherm_other,"Second parameters (alpha, ...) defining isotherm c_s = omega * (alpha*c_a)/(1- alpha*c_a).","1.0");
66  ADD_FIELD(init_conc_solid, "Initial solid concentration of substances."
67  " Vector, one value for every substance.", "0");
68 
69  rock_density.units("");
70  init_conc_solid.units("M/L^3");
71 
72  input_data_set_ += *this;
73 
74  // porosity field is set from governing equation (transport) later
75  // hence we do not add it to the input_data_set_
76  *this += porosity.name("porosity").units("1").just_copy();
77 
78  output_fields += *this;
79  output_fields += conc_solid.name(output_field_name).units("M/L^3");
80 }
81 
82 
84  : ReactionTerm(init_mesh, in_rec),
85  data_(nullptr)
86 {
87 }
88 
89 
91 {
92  if(reaction != nullptr) delete reaction;
93 
94  if (data_ != nullptr) delete data_;
95 
96 // if(!output_rec.is_empty())
97  {
98  VecDestroy(vconc_solid);
99  VecDestroy(vconc_solid_out);
100  }
101 
102  for (unsigned int sbi = 0; sbi < names_.size(); sbi++)
103  {
104  //no mpi vectors
105  xfree(conc_solid[sbi]);
106  }
107  xfree(conc_solid);
108 }
109 
110 
111 
112 
114 {
115  Input::Array substances_array = in_rec.val<Input::Array>("substances");
116  unsigned int k, idx, i_spec = 0;
117 
118  for(Input::Iterator<string> spec_iter = substances_array.begin<string>(); spec_iter != substances_array.end(); ++spec_iter, i_spec++)
119  {
120  //finding name in the global array of names
121  for(k = 0; k < names.size(); k++)
122  {
123  if (*spec_iter == names[k])
124  {
125  idx = k;
126  break;
127  }
128  }
129 
130  if ((idx < names.size()) && (idx >= 0))
131  {
132  substance_id[i_spec] = idx; //mapping - if not found, it creates new map
133  }
134  else xprintf(UsrErr,"Wrong name of %d-th reaction specie - not found in global set of transported substances.\n", i_spec);
135  }
136  n_substances_ = substance_id.size();
137 }
138 
140 {
145 
146  // read fields from input file
147  data_->input_data_set_.set_input_list(in_rec.val<Input::Array>("input_fields"));
148 
149  data_->set_mesh(*mesh_);
150  data_->set_limit_side(LimitSide::right);
151 
152  // Common data for all the isotherms loaded bellow
153  solvent_density = in_rec.val<double>("solvent_density");
154 
155  Input::Array molar_mass_array = in_rec.val<Input::Array>("molar_mass");
156 
157  if (molar_mass_array.size() == molar_masses.size() ) molar_mass_array.copy_to( molar_masses );
158  else xprintf(UsrErr,"Number of molar masses %d has to match number of adsorbing species %d.\n", molar_mass_array.size(), molar_masses.size());
159 
160  Input::Iterator<Input::Array> solub_iter = in_rec.find<Input::Array>("solubility");
161  if( solub_iter )
162  {
163  solub_iter->copy_to(solubility_vec_);
164  if (solubility_vec_.size() != n_substances_)
165  {
166  xprintf(UsrErr,"Number of given solubility limits %d has to match number of adsorbing species %d.\n", solubility_vec_.size(), n_substances_);
167  }
168  }else{
169  // fill solubility_vec_ with zeros or resize it at least
171  }
172 
173  Input::Iterator<Input::Array> interp_table_limits = in_rec.find<Input::Array>("table_limits");
174  if( interp_table_limits )
175  {
176  interp_table_limits->copy_to(table_limit_);
177  if (table_limit_.size() != n_substances_)
178  {
179  xprintf(UsrErr,"Number of given table limits %d has to match number of adsorbing species %d.\n", table_limit_.size(), n_substances_);
180  }/**/
181  }else{
182  // fill table_limit_ with zeros or resize it at least
183  table_limit_.resize(n_substances_);
184  }
185 }
186 
188 {
189  ASSERT(distribution != nullptr, "Distribution has not been set yet.\n");
190  ASSERT(time_ != nullptr, "Time governor has not been set yet.\n");
191 
192 
193 
194  //DBGMSG("SorptionBase constructor.\n");
196 
197  // for(unsigned int s=0; s<names_.size(); s++)
198  // cout << s << " " << names_[s] << endl;
199  //
200  // for(unsigned int s=0; s<n_substances_; s++)
201  // cout << s << " " << substance_id[s] << " " << names_[substance_id[s]] << endl;
202 
204  nr_of_points = input_record_.val<int>("substeps");
205 
206  // Input::Iterator<Input::Record> out_rec = in_rec.find<Input::Record>("output");
207  //output_rec = in_rec.find<Input::Record>("output");
208  // if(out_rec) output_rec = *out_rec;
209  output_array = input_record_.val<Input::Array>("output_fields");
210 
211  //Simple vectors holding common informations.
212  molar_masses.resize( n_substances_ );
213 
214  //isotherms array resized bellow
215  isotherms.resize(nr_of_regions);
216  for(int i_reg = 0; i_reg < nr_of_regions; i_reg++)
217  for(int i_spec = 0; i_spec < n_substances_; i_spec++)
218  {
219  Isotherm iso_mob;
220  isotherms[i_reg].push_back(iso_mob);
221  }
222 
223 
225  data_->set_time(*time_);
226  make_tables();
227 
228  //allocating new array for sorbed concentrations
229  conc_solid = (double**) xmalloc(names_.size() * sizeof(double*));//new double * [n_substances_];
230  conc_solid_out = (double**) xmalloc(names_.size() * sizeof(double*));
231  for (unsigned int sbi = 0; sbi < names_.size(); sbi++)
232  {
233  conc_solid[sbi] = (double*) xmalloc(distribution->lsize() * sizeof(double));//new double[ nr_of_local_elm ];
234  conc_solid_out[sbi] = (double*) xmalloc(distribution->size() * sizeof(double));
235  //zero initialization of solid concentration for all substances
236  for(unsigned int i=0; i < distribution->lsize(); i++)
237  conc_solid[sbi][i] = 0;
238  }
239 
241 
242 
243  //setting initial condition for solid concentrations
244  for (unsigned int loc_el = 0; loc_el < distribution->lsize(); loc_el++)
245  {
246  unsigned int index = el_4_loc[loc_el];
247  ElementAccessor<3> ele_acc = mesh_->element_accessor(index);
248  arma::vec value = data_->init_conc_solid.value(ele_acc.centre(),
249  ele_acc);
250 
251  //setting initial solid concentration for substances involved in adsorption
252  for (int sbi = 0; sbi < n_substances_; sbi++)
253  {
254  int subst_id = substance_id[sbi];
255  conc_solid[subst_id][loc_el] = value(sbi);
256  }
257  }
258 
259  //initialization of output
263  for (int sbi=0; sbi<names_.size(); sbi++)
264  {
265  // create shared pointer to a FieldElementwise and push this Field to output_field on all regions
266  std::shared_ptr<FieldElementwise<3, FieldValue<3>::Scalar> > output_field_ptr(
268  data_->conc_solid[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
269  }
270  data_->output_fields.set_limit_side(LimitSide::right);
272 
273 
274  //DBGMSG("Going to write initial condition.\n");
275  // write initial condition
279 
280  // creating reaction from input and setting their parameters
282 
283  if(reaction != nullptr)
284  {
289  }
290 }
291 
293 {
294  //DBGMSG("SorptionBase init_from_input\n");
295  Input::Iterator<Input::AbstractRecord> reactions_it = in_rec.find<Input::AbstractRecord>("reaction");
296  if ( reactions_it )
297  {
298  if (reactions_it->type() == Linear_reaction::input_type ) {
299  reaction = new Linear_reaction(*mesh_, *reactions_it);
300 
301  } else
302  if (reactions_it->type() == Pade_approximant::input_type) {
303  reaction = new Pade_approximant(*mesh_, *reactions_it);
304 
305  } else
306  if (reactions_it->type() == SorptionBase::input_type ) {
307  xprintf(UsrErr, "Sorption model cannot have another descendant sorption model.\n");
308  } else
309  if (reactions_it->type() == DualPorosity::input_type ) {
310  xprintf(UsrErr, "Sorption model cannot have descendant dual porosity model.\n");
311  } else
312  if (reactions_it->type() == Semchem_interface::input_type )
313  {
314  xprintf(UsrErr, "Semchem chemistry model is not supported at current time.\n");
315  } else
316  {
317  xprintf(UsrErr, "Unknown reactions type in Sorption model.\n");
318  }
319  } else
320  {
321  reaction = nullptr;
322  }
323 }
324 
326 {
327  //DBGMSG("Sorption - update_solution\n");
328  data_->set_time(*time_); // set to the last computed time
329 
330  // if parameters changed during last time step, reinit isotherms and eventualy
331  // update interpolation tables in the case of constant rock matrix parameters
332  if(data_->changed())
333  make_tables();
334 
335 
336  START_TIMER("Sorption");
337  for (int loc_el = 0; loc_el < distribution->lsize(); loc_el++)
338  {
340  }
341  END_TIMER("Sorption");
342 
343  if(reaction != nullptr) reaction->update_solution();
344 }
345 
347 {
348  ElementAccessor<3> elm;
349  BOOST_FOREACH(const Region &reg_iter, this->mesh_->region_db().get_region_set("BULK") )
350  {
351  int reg_idx = reg_iter.bulk_idx();
352 
353  if(data_->is_constant(reg_iter))
354  {
355  ElementAccessor<3> elm(this->mesh_, reg_iter); // constant element accessor
356  isotherm_reinit(isotherms[reg_idx],elm);
357  //xprintf(MsgDbg,"parameters are constant\n");
358  for(int i_subst = 0; i_subst < n_substances_; i_subst++)
359  {
360  isotherms[reg_idx][i_subst].make_table(nr_of_points);
361  }
362  }
363  }
364 }
365 
366 double **SorptionBase::compute_reaction(double **concentrations, int loc_el) // Sorption simulations are realized just for one element.
367 {
368  //DBGMSG("compute_reaction\n");
369  ElementFullIter elem = mesh_->element(el_4_loc[loc_el]);
370  double porosity;
371  double rock_density;
372  Region region = elem->region();
373  int reg_id_nr = region.bulk_idx();
374  int variabl_int = 0;
375  int i_subst, subst_id;
376 
377  std::vector<Isotherm> & isotherms_vec = isotherms[reg_id_nr];
378 
379  //if(reg_id_nr != 0) cout << "region id is " << reg_id_nr << endl;
380 
381  // Constant value of rock density and mobile porosity over the whole region => interpolation_table is precomputed
382  if (isotherms_vec[0].is_precomputed())
383  {
384  for(i_subst = 0; i_subst < n_substances_; i_subst++)
385  {
386  subst_id = substance_id[i_subst];
387  //DBGMSG("on s_%d precomputed %d\n",subst_id, isotherms_vec[i_subst].is_precomputed());
388 
389  Isotherm & isotherm = this->isotherms[reg_id_nr][i_subst];
390  isotherm.interpolate((concentration_matrix_[subst_id][loc_el]), conc_solid[subst_id][loc_el]);
391  }
392  }
393  else
394  {
395  isotherm_reinit(isotherms_vec, elem->element_accessor());
396 
397  for(i_subst = 0; i_subst < n_substances_; i_subst++)
398  {
399  subst_id = substance_id[i_subst];
400  Isotherm & isotherm = this->isotherms[reg_id_nr][i_subst];
401  isotherm.compute((concentration_matrix_[subst_id][loc_el]), conc_solid[subst_id][loc_el]);
402  }
403  }
404 
405  return concentrations;
406 }
407 
409 {
410  data_->set_field(data_->porosity.name(),por_m);
411 }
412 
413 
415 {
416  DBGMSG("Not implemented.\n");
417  xprintf(Msg, "\nSorption parameters are defined as follows:\n");
418 }
419 
420 /*
421 void SorptionBase::set_concentration_vector(Vec &vc)
422 {
423  DBGMSG("Not implemented.\n");
424 }
425 */
426 
427 /**************************************** OUTPUT ***************************************************/
428 
430 {
431  int sbi, n_subst, ierr, rank, np; //, i, j, ph;
432  n_subst = names_.size();
433 
434  vconc_solid = (Vec*) xmalloc(n_subst * (sizeof(Vec)));
435  vconc_solid_out = (Vec*) xmalloc(n_subst * (sizeof(Vec))); // extend to all
436 
437  for (sbi = 0; sbi < n_subst; sbi++) {
438  ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,1, distribution->lsize(), mesh_->n_elements(), conc_solid[sbi],
439  &vconc_solid[sbi]);
440  VecZeroEntries(vconc_solid[sbi]);
441 
442  ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1, mesh_->n_elements(), conc_solid_out[sbi], &vconc_solid_out[sbi]);
443  VecZeroEntries(vconc_solid_out[sbi]);
444  }
445 }
446 
447 
449 {
450  unsigned int sbi/*, rank, np*/;
451  IS is;
452  VecScatter vconc_out_scatter;
453  //PetscViewer inviewer;
454 
455  ISCreateGeneral(PETSC_COMM_SELF, mesh_->n_elements(), row_4_el, PETSC_COPY_VALUES, &is); //WithArray
456  VecScatterCreate(vconc_solid[0], is, vconc_solid_out[0], PETSC_NULL, &vconc_out_scatter);
457  for (sbi = 0; sbi < names_.size(); sbi++) {
458  VecScatterBegin(vconc_out_scatter, vconc_solid[sbi], vconc_solid_out[sbi], INSERT_VALUES, SCATTER_FORWARD);
459  VecScatterEnd(vconc_out_scatter, vconc_solid[sbi], vconc_solid_out[sbi], INSERT_VALUES, SCATTER_FORWARD);
460  }
461  //VecView(transport->vconc[0],PETSC_VIEWER_STDOUT_WORLD);
462  //VecView(transport->vconc_out[0],PETSC_VIEWER_STDOUT_WORLD);
463  VecScatterDestroy(&(vconc_out_scatter));
464  ISDestroy(&(is));
465 }
466 
467 
469 {
470  //DBGMSG("Sorption output\n");
472 
473  int ierr, rank;
474  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
475  if (rank == 0)
476  {
477  // Register fresh output data
480  }
481 
482  //it can call only linear reaction which has no output at the moment
483  //if(reaction) reaction->output_data();
484 
485  //for synchronization when measuring time by Profiler
487 }