Flow123d  last_with_con_2.0.0-4-g42e6930
dual_porosity.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 dual_porosity.cc
15  * @brief
16  */
17 
18 #include <iostream>
19 #include <stdlib.h>
20 #include <math.h>
21 
24 #include "system/system.hh"
25 #include "system/sys_profiler.hh"
26 
27 #include "la/distribution.hh"
28 #include "mesh/mesh.h"
29 #include "mesh/elements.h"
30 #include "mesh/region.hh"
32 
33 #include "reaction/sorption.hh"
37 #include "input/factory.hh"
38 
39 FLOW123D_FORCE_LINK_IN_CHILD(dualPorosity)
40 
41 
42 using namespace Input::Type;
43 
44 
45 
46 const Record & DualPorosity::get_input_type() {
47  return Record("DualPorosity",
48  "Dual porosity model in transport problems.\n"
49  "Provides computing the concentration of substances in mobile and immobile zone.\n"
50  )
52  .declare_key("input_fields", Array(DualPorosity::EqData().make_field_descriptor_type("DualPorosity")), Default::obligatory(),
53  "Containes region specific data necessary to construct dual porosity model.")
54  .declare_key("scheme_tolerance", Double(0.0), Default("1e-3"),
55  "Tolerance according to which the explicit Euler scheme is used or not."
56  "Set 0.0 to use analytic formula only (can be slower).")
57 
58  .declare_key("reaction_mobile", ReactionTerm::get_input_type(), Default::optional(), "Reaction model in mobile zone.")
59  .declare_key("reaction_immobile", ReactionTerm::get_input_type(), Default::optional(), "Reaction model in immobile zone.")
60  .declare_key("output",
61  EqData().output_fields.make_output_type("DualPorosity", ""),
62  IT::Default("{ \"fields\": [ \"conc_immobile\" ] }"),
63  "Setting of the fields output.")
64  .close();
65 }
66 
67 const int DualPorosity::registrar =
68  Input::register_class< DualPorosity, Mesh &, Input::Record >("DualPorosity") +
70 
72 {
74  .name("diffusion_rate_immobile")
75  .description("Diffusion coefficient of non-equilibrium linear exchange between mobile and immobile zone.")
76  .input_default("0")
77  .units( UnitSI().s(-1) );
78 
79  *this += porosity_immobile
80  .name("porosity_immobile")
81  .description("Porosity of the immobile zone.")
82  .input_default("0")
84 
85  *this += init_conc_immobile
86  .name("init_conc_immobile")
87  .description("Initial concentration of substances in the immobile zone.")
88  .units( UnitSI().kg().m(-3) );
89 
90  //creating field for porosity that is set later from the governing equation (transport)
91  *this +=porosity
92  .name("porosity")
95 
96  *this += conc_immobile
97  .name("conc_immobile")
98  .units( UnitSI().kg().m(-3) )
100 
101  output_fields += *this;
102 
103 }
104 
106  : ReactionTerm(init_mesh, in_rec)
107 {
108  //set pointer to equation data fieldset
109  this->eq_data_ = &data_;
110 
111  //reads input and creates possibly other reaction terms
112  make_reactions();
113  //read value from input
114  scheme_tolerance_ = input_record_.val<double>("scheme_tolerance");
115 }
116 
118 {
119  VecScatterDestroy(&(vconc_out_scatter));
120 
121 
122  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
123  {
124 
125  //no mpi vectors
126  VecDestroy(&(vconc_immobile[sbi]));
127  delete [] conc_immobile[sbi];
128  }
129 
130  delete [] vconc_immobile;
131  delete [] conc_immobile;
132 }
133 
134 
137  if ( reactions_it )
138  {
139  // TODO: allowed instances in this case are only
140  // FirstOrderReaction, RadioactiveDecay and SorptionMob
141  reaction_mobile = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
142  } else
143  {
144  reaction_mobile = nullptr;
145  }
146 
147  reactions_it = input_record_.find<Input::AbstractRecord>("reaction_immobile");
148  if ( reactions_it )
149  {
150  // TODO: allowed instances in this case are only
151  // FirstOrderReaction, RadioactiveDecay and SorptionImmob
152  reaction_immobile = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
153  } else
154  {
155  reaction_immobile = nullptr;
156  }
157 
158 }
159 
161 {
162  OLD_ASSERT(distribution_ != nullptr, "Distribution has not been set yet.\n");
163  OLD_ASSERT(time_ != nullptr, "Time governor has not been set yet.\n");
164  OLD_ASSERT(output_stream_,"Null output stream.");
166 
167  //allocating memory for immobile concentration matrix
168  conc_immobile = new double* [substances_.size()];
169  conc_immobile_out.clear();
170  conc_immobile_out.resize( substances_.size() );
171  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
172  {
173  conc_immobile[sbi] = new double [distribution_->lsize()];
174  conc_immobile_out[sbi].resize( distribution_->size() );
175  }
177 
179 
180  if(reaction_mobile)
181  {
182  reaction_mobile->substances(substances_)
183  .output_stream(output_stream_)
184  .concentration_matrix(concentration_matrix_, distribution_, el_4_loc_, row_4_el_)
185  .set_time_governor(*time_);
186  reaction_mobile->initialize();
187  }
188 
190  {
191  reaction_immobile->substances(substances_)
192  .output_stream(output_stream_)
193  .concentration_matrix(conc_immobile, distribution_, el_4_loc_, row_4_el_)
194  .set_time_governor(*time_);
195  reaction_immobile->initialize();
196  }
197 
198 }
199 
201 {
203  //setting fields that are set from input file
206 
207  //setting fields in data
208  data_.set_mesh(*mesh_);
209 
210  //initialization of output
211  //output_array = input_record_.val<Input::Array>("output_fields");
216  for (unsigned int sbi=0; sbi<substances_.size(); sbi++)
217  {
218  // create shared pointer to a FieldElementwise and push this Field to output_field on all regions
219  auto output_field_ptr = conc_immobile_out[sbi].create_field<3, FieldValue<3>::Scalar>(substances_.size());
220  data_.conc_immobile[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr, 0);
221  }
222  //output_stream_->add_admissible_field_names(output_array);
224 }
225 
226 
228 {
229  OLD_ASSERT(distribution_ != nullptr, "Distribution has not been set yet.\n");
230  OLD_ASSERT(time_ != nullptr, "Time governor has not been set yet.\n");
231  OLD_ASSERT(output_stream_,"Null output stream.");
233 
234  //coupling - passing fields
235  if(reaction_mobile)
236  if (typeid(*reaction_mobile) == typeid(SorptionMob))
237  {
238  reaction_mobile->data().set_field("porosity", data_["porosity"]);
239  reaction_mobile->data().set_field("porosity_immobile", data_["porosity_immobile"]);
240  }
242  if (typeid(*reaction_immobile) == typeid(SorptionImmob))
243  {
244  reaction_immobile->data().set_field("porosity", data_["porosity"]);
245  reaction_immobile->data().set_field("porosity_immobile", data_["porosity_immobile"]);
246  }
247 
250 
251  // write initial condition
252  //output_vector_gather();
253  //data_.output_fields.set_time(time_->step(0), LimitSide::right);
254  //data_.output_fields.output(output_stream_);
255 
256  output_data();
257 
258  if(reaction_mobile)
259  reaction_mobile->zero_time_step();
260 
262  reaction_immobile->zero_time_step();
263 
264 }
265 
267 {
268  //setting initial condition for immobile concentration matrix
269  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
270  {
271  unsigned int index = el_4_loc_[loc_el];
272  ElementAccessor<3> ele_acc = mesh_->element_accessor(index);
273  arma::vec value = data_.init_conc_immobile.value(ele_acc.centre(), ele_acc);
274 
275  for (unsigned int sbi=0; sbi < substances_.size(); sbi++)
276  {
277  conc_immobile[sbi][loc_el] = value(sbi);
278  }
279  }
280 }
281 
283 {
285 
286  START_TIMER("dual_por_exchange_step");
287  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
288  {
290  }
291  END_TIMER("dual_por_exchange_step");
292 
293  if(reaction_mobile) reaction_mobile->update_solution();
294  if(reaction_immobile) reaction_immobile->update_solution();
295 }
296 
297 
298 double **DualPorosity::compute_reaction(double **concentrations, int loc_el)
299 {
300  unsigned int sbi;
301  double conc_average, // weighted (by porosity) average of concentration
302  conc_mob, conc_immob, // new mobile and immobile concentration
303  previous_conc_mob, previous_conc_immob, // mobile and immobile concentration in previous time step
304  conc_max, //difference between concentration and average concentration
305  por_mob, por_immob; // mobile and immobile porosity
306 
307  // get data from fields
308  ElementFullIter ele = mesh_->element(el_4_loc_[loc_el]);
309  por_mob = data_.porosity.value(ele->centre(),ele->element_accessor());
310  por_immob = data_.porosity_immobile.value(ele->centre(),ele->element_accessor());
311  arma::Col<double> diff_vec = data_.diffusion_rate_immobile.value(ele->centre(), ele->element_accessor());
312 
313  // if porosity_immobile == 0 then mobile concentration stays the same
314  // and immobile concentration cannot change
315  if (por_immob == 0.0) return conc_immobile;
316 
317  double exponent,
318  temp_exponent = (por_mob + por_immob) / (por_mob * por_immob) * time_->dt();
319 
320  for (sbi = 0; sbi < substances_.size(); sbi++) //over all substances
321  {
322  exponent = diff_vec[sbi] * temp_exponent;
323  //previous values
324  previous_conc_mob = concentration_matrix_[sbi][loc_el];
325  previous_conc_immob = conc_immobile[sbi][loc_el];
326 
327  // ---compute average concentration------------------------------------------
328  conc_average = ((por_mob * previous_conc_mob) + (por_immob * previous_conc_immob))
329  / (por_mob + por_immob);
330 
331  conc_max = std::max(previous_conc_mob-conc_average, previous_conc_immob-conc_average);
332 
333  // the following 2 conditions guarantee:
334  // 1) stability of forward Euler's method
335  // 2) the error of forward Euler's method will not be large
336  if(time_->dt() <= por_mob*por_immob/(max(diff_vec)*(por_mob+por_immob)) &&
337  conc_max <= (2*scheme_tolerance_/(exponent*exponent)*conc_average)) // forward euler
338  {
339  double temp = diff_vec[sbi]*(previous_conc_immob - previous_conc_mob) * time_->dt();
340  // ---compute concentration in mobile area
341  conc_mob = temp / por_mob + previous_conc_mob;
342 
343  // ---compute concentration in immobile area
344  conc_immob = -temp / por_immob + previous_conc_immob;
345  }
346  else //analytic solution
347  {
348  double temp = exp(-exponent);
349  // ---compute concentration in mobile area
350  conc_mob = (previous_conc_mob - conc_average) * temp + conc_average;
351 
352  // ---compute concentration in immobile area
353  conc_immob = (previous_conc_immob - conc_average) * temp + conc_average;
354  }
355 
356  concentration_matrix_[sbi][loc_el] = conc_mob;
357  conc_immobile[sbi][loc_el] = conc_immob;
358  }
359 
360  return conc_immobile;
361 }
362 
363 
365 {
366  int sbi, n_subst;
367  n_subst = substances_.size();
368 
369  vconc_immobile = new Vec [n_subst];
370 
371 
372  for (sbi = 0; sbi < n_subst; sbi++) {
373  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, distribution_->lsize(), mesh_->n_elements(), conc_immobile[sbi],
374  &vconc_immobile[sbi]);
375  VecZeroEntries(vconc_immobile[sbi]);
376 
377  // if(rank == 0)
378  VecZeroEntries(conc_immobile_out[sbi].get_data_petsc());
379  }
380 
381  // create output vector scatter
382  IS is;
383  ISCreateGeneral(PETSC_COMM_SELF, mesh_->n_elements(), row_4_el_, PETSC_COPY_VALUES, &is); //WithArray
384  VecScatterCreate(vconc_immobile[0], is, conc_immobile_out[0].get_data_petsc(), PETSC_NULL, &vconc_out_scatter);
385  ISDestroy(&(is));
386 }
387 
388 
390 {
391  unsigned int sbi;
392 
393  for (sbi = 0; sbi < substances_.size(); sbi++) {
394  VecScatterBegin(vconc_out_scatter, vconc_immobile[sbi], conc_immobile_out[sbi].get_data_petsc(), INSERT_VALUES, SCATTER_FORWARD);
395  VecScatterEnd(vconc_out_scatter, vconc_immobile[sbi], conc_immobile_out[sbi].get_data_petsc(), INSERT_VALUES, SCATTER_FORWARD);
396  }
397 }
398 
399 
401 {
405  }
406 
407  // Register fresh output data
409 
410  if (time_->tlevel() !=0) {
411  // zero_time_step call zero_time_Step of subreactions which performs its own output
412  if (reaction_mobile) reaction_mobile->output_data();
413  if (reaction_immobile) reaction_immobile->output_data();
414  }
415 }
416 
417 
418 bool DualPorosity::evaluate_time_constraint(double &time_constraint)
419 {
420  bool cfl_changed = false;
421  if (reaction_mobile)
422  {
423  if (reaction_mobile->evaluate_time_constraint(time_constraint))
424  cfl_changed = true;
425  }
426  if (reaction_immobile)
427  {
428  double cfl_immobile;
429  if (reaction_immobile->evaluate_time_constraint(cfl_immobile))
430  {
431  time_constraint = std::min(time_constraint, cfl_immobile);
432  cfl_changed = true;
433  }
434  }
435 
436  return cfl_changed;
437 }
TimeGovernor & time()
Definition: equation.hh:146
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:203
unsigned int size() const
get global size
void set_input_list(Input::Array input_list)
Definition: field_set.hh:180
Sorption model in immobile zone in case dual porosity is considered.
Definition: sorption.hh:123
FieldSet * eq_data_
Definition: equation.hh:230
void set_initial_condition()
Sets initial condition from input.
VecScatter vconc_out_scatter
Output vector scatter.
std::shared_ptr< ReactionTerm > reaction_mobile
Reaction running in mobile zone.
Accessor to input data conforming to declared Array.
Definition: accessors.hh:552
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:582
int tlevel() const
void output_data(void) override
Main output routine.
double ** concentration_matrix_
const std::vector< std::string > & names()
Definition: substance.hh:85
MultiField< 3, FieldValue< 3 >::Scalar > conc_immobile
Calculated concentrations in the immobile zone.
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:50
static Input::Type::Abstract & get_input_type()
void set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.hh:195
void output(TimeStep step)
SubstanceList substances_
Names belonging to substances.
EquationOutput output_fields
Fields indended for output, i.e. all input fields plus those representing solution.
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
void zero_time_step() override
Definition: mesh.h:95
Field< 3, FieldValue< 3 >::Scalar > porosity_immobile
Immobile porosity field.
Iterator< Ret > find(const string &key) const
int * row_4_el_
Indices of rows belonging to elements.
void initialize(std::shared_ptr< OutputTime > stream, Input::Record in_rec, const TimeGovernor &tg)
~DualPorosity(void)
Destructor.
std::vector< VectorSeqDouble > conc_immobile_out
concentration array output for immobile phase (gathered - sequential)
const RegionDB & region_db() const
Definition: mesh.h:152
ElementAccessor< 3 > element_accessor(unsigned int idx, bool boundary=false)
Definition: mesh.cc:700
void allocate_output_mpi(void)
Allocates petsc vectors, prepares them for output and creates output vector scatter.
const TimeStep & step(int index=-1) const
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
unsigned int size() const
Definition: substance.hh:87
Class for declaration of inputs sequences.
Definition: type_base.hh:321
std::shared_ptr< OutputTime > output_stream_
Pointer to a transport output stream.
void update_solution(void) override
ReactionTerm(Mesh &init_mesh, Input::Record in_rec)
Constructor.
Class ReactionTerm is an abstract class representing reaction term in transport.
Class Dual_por_exchange implements the model of dual porosity.
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:195
double scheme_tolerance_
Dual porosity computational scheme tolerance.
#define OLD_ASSERT(...)
Definition: global_defs.h:131
void setup_components()
std::shared_ptr< ReactionTerm > reaction_immobile
Reaction running in immobile zone.
void make_reactions()
Resolves construction of following reactions.
unsigned int n_elements() const
Definition: mesh.h:138
void initialize() override
Prepares the object to usage.
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:518
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
static constexpr Mask input_copy
Definition: field_flag.hh:44
FieldCommon & input_default(const string &input_default)
Accessor to the data with type Type::Record.
Definition: accessors.hh:277
const Ret val(const string &key) const
#define START_TIMER(tag)
Starts a timer with specified tag.
#define OLD_ASSERT_LESS(a, b)
Definition: global_defs.h:134
bool evaluate_time_constraint(double &time_constraint) override
Computes a constraint for time step.
void initialize_fields()
Initializes field sets.
Mesh * mesh_
Definition: equation.hh:221
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:342
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:468
EqData()
Collect all fields.
Class representing dual porosity model in transport.
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:444
double ** conc_immobile
Support classes for parallel programing.
Sorption model in mobile zone in case dual porosity is considered.
Definition: sorption.hh:95
void output_vector_gather(void) override
Gathers all the parallel vectors to enable them to be output.
FieldCommon & description(const string &description)
Definition: field_common.hh:93
virtual MultiFieldValue::return_type value(const Point &p, const ElementAccessor< spacedim > &elm) const
Input::Record input_record_
Definition: equation.hh:223
void set_components(const std::vector< string > &names)
Definition: field_set.hh:167
int * el_4_loc_
Indices of elements belonging to local dofs.
double dt() const
double ** compute_reaction(double **concentrations, int loc_el) override
FieldCommon & name(const string &name)
Definition: field_common.hh:86
#define END_TIMER(tag)
Ends a timer with specified tag.
MultiField< 3, FieldValue< 3 >::Scalar > diffusion_rate_immobile
Mass transfer coefficients between mobile and immobile pores.
arma::vec::fixed< spacedim > centre() const
Definition: accessors.hh:91
static const int registrar
Registrar of class to factory.
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:173
MultiField< 3, FieldValue< 3 >::Scalar > init_conc_immobile
Initial concentrations in the immobile zone.
Record type proxy class.
Definition: type_record.hh:171
bool is_field_output_time(const FieldCommon &field, TimeStep step) const
Field< 3, FieldValue< 3 >::Scalar > porosity
Porosity field.
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Class for representation SI units of Fields.
Definition: unit_si.hh:40
static const Input::Type::Record & get_input_type()
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:53
DualPorosity data.
Vec * vconc_immobile
PETSC concentration vector for immobile phase (parallel).
This file contains classes representing sorption model. Sorption model can be computed both in case t...
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
Definition: global_defs.h:180
TimeGovernor * time_
Definition: equation.hh:222
FieldSet input_data_set_
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:221
unsigned int lsize(int proc) const
get local size