Flow123d  release_3.0.0-684-g928e266
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"
31 #include "mesh/accessors.hh"
32 #include "fields/field_fe.hh"
34 
35 #include "reaction/sorption.hh"
38 #include "input/factory.hh"
39 
40 FLOW123D_FORCE_LINK_IN_CHILD(dualPorosity)
41 
42 
43 using namespace Input::Type;
44 
45 
46 
47 const Record & DualPorosity::get_input_type() {
48  return Record("DualPorosity",
49  "Dual porosity model in transport problems.\n"
50  "Provides computing the concentration of substances in mobile and immobile zone.\n"
51  )
53  .declare_key("input_fields", Array(DualPorosity::EqData().make_field_descriptor_type("DualPorosity")), Default::obligatory(),
54  "Containes region specific data necessary to construct dual porosity model.")
55  .declare_key("scheme_tolerance", Double(0.0), Default("1e-3"),
56  "Tolerance according to which the explicit Euler scheme is used or not."
57  "Set 0.0 to use analytic formula only (can be slower).")
58 
59  .declare_key("reaction_mobile", ReactionTerm::it_abstract_mobile_term(), Default::optional(), "Reaction model in mobile zone.")
60  .declare_key("reaction_immobile", ReactionTerm::it_abstract_immobile_term(), Default::optional(), "Reaction model in immobile zone.")
61  .declare_key("output",
62  EqData().output_fields.make_output_type("DualPorosity", ""),
63  IT::Default("{ \"fields\": [ \"conc_immobile\" ] }"),
64  "Setting of the fields output.")
65  .close();
66 }
67 
68 const int DualPorosity::registrar =
69  Input::register_class< DualPorosity, Mesh &, Input::Record >("DualPorosity") +
71 
73 {
75  .name("diffusion_rate_immobile")
76  .description("Diffusion coefficient of non-equilibrium linear exchange between mobile and immobile zone.")
77  .input_default("0")
78  .units( UnitSI().s(-1) );
79 
80  *this += porosity_immobile
81  .name("porosity_immobile")
82  .description("Porosity of the immobile zone.")
83  .input_default("0")
85  .set_limits(0.0);
86 
87  *this += init_conc_immobile
88  .name("init_conc_immobile")
89  .description("Initial concentration of substances in the immobile zone.")
90  .units( UnitSI().kg().m(-3) );
91 
92  //creating field for porosity that is set later from the governing equation (transport)
93  *this +=porosity
94  .name("porosity")
95  .description("Concentration solution in the mobile phase.")
98  .set_limits(0.0);
99 
100  *this += conc_immobile
101  .name("conc_immobile")
102  .units( UnitSI().kg().m(-3) )
104 
105  output_fields += *this;
106 
107 }
108 
110  : ReactionTerm(init_mesh, in_rec)
111 {
112  //set pointer to equation data fieldset
113  this->eq_data_ = &data_;
114 
115  //reads input and creates possibly other reaction terms
116  make_reactions();
117  //read value from input
118  scheme_tolerance_ = input_record_.val<double>("scheme_tolerance");
119 }
120 
122 {
123  chkerr(VecScatterDestroy(&(vconc_out_scatter)));
124 
125 
126  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
127  {
128 
129  //no mpi vectors
130  chkerr(VecDestroy( &vconc_immobile[sbi] ));
131  delete [] conc_immobile[sbi];
132  }
133 
134  delete [] vconc_immobile;
135  delete [] conc_immobile;
136 }
137 
138 
141  if ( reactions_it )
142  {
143  // TODO: allowed instances in this case are only
144  // FirstOrderReaction, RadioactiveDecay and SorptionMob
145  reaction_mobile = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
146  } else
147  {
148  reaction_mobile = nullptr;
149  }
150 
151  reactions_it = input_record_.find<Input::AbstractRecord>("reaction_immobile");
152  if ( reactions_it )
153  {
154  // TODO: allowed instances in this case are only
155  // FirstOrderReaction, RadioactiveDecay and SorptionImmob
156  reaction_immobile = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(*mesh_, *reactions_it);
157  } else
158  {
159  reaction_immobile = nullptr;
160  }
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  //allocating memory for immobile concentration matrix
172  conc_immobile = new double* [substances_.size()];
173  conc_immobile_out.clear();
174  conc_immobile_out.resize( substances_.size() );
175  output_field_ptr.clear();
176  output_field_ptr.resize( substances_.size() );
177  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++)
178  {
179  conc_immobile[sbi] = new double [distribution_->lsize()];
180  conc_immobile_out[sbi].resize( distribution_->size() );
181  }
183 
185 
186  if(reaction_mobile)
187  {
188  reaction_mobile->substances(substances_)
189  .output_stream(output_stream_)
190  .concentration_matrix(concentration_matrix_, distribution_, el_4_loc_, row_4_el_)
191  .set_time_governor(*time_);
192  reaction_mobile->initialize();
193  }
194 
196  {
197  reaction_immobile->substances(substances_)
198  .output_stream(output_stream_)
199  .concentration_matrix(conc_immobile, distribution_, el_4_loc_, row_4_el_)
200  .set_time_governor(*time_);
201  reaction_immobile->initialize();
202  }
203 
204 }
205 
207 {
209  //setting fields that are set from input file
212 
213  //setting fields in data
214  data_.set_mesh(*mesh_);
215 
216  //initialization of output
221  for (unsigned int sbi=0; sbi<substances_.size(); sbi++)
222  {
223  // create shared pointer to a FieldFE and push this Field to output_field on all regions
224  output_field_ptr[sbi] = create_field<3, FieldValue<3>::Scalar>(conc_immobile_out[sbi], *mesh_, 1);
225  data_.conc_immobile[sbi].set_field(mesh_->region_db().get_region_set("ALL"), output_field_ptr[sbi], 0);
226  }
228 }
229 
230 
232 {
233  OLD_ASSERT(distribution_ != nullptr, "Distribution has not been set yet.\n");
234  OLD_ASSERT(time_ != nullptr, "Time governor has not been set yet.\n");
235  OLD_ASSERT(output_stream_,"Null output stream.");
237 
238  //coupling - passing fields
239  if(reaction_mobile)
240  if (typeid(*reaction_mobile) == typeid(SorptionMob))
241  {
242  reaction_mobile->data().set_field("porosity", data_["porosity"]);
243  reaction_mobile->data().set_field("porosity_immobile", data_["porosity_immobile"]);
244  }
246  if (typeid(*reaction_immobile) == typeid(SorptionImmob))
247  {
248  reaction_immobile->data().set_field("porosity", data_["porosity"]);
249  reaction_immobile->data().set_field("porosity_immobile", data_["porosity_immobile"]);
250  }
251 
253  std::stringstream ss; // print warning message with table of uninitialized fields
254  if ( FieldCommon::print_message_table(ss, "dual porosity") ) {
255  WarningOut() << ss.str();
256  }
258 
259  output_data();
260 
261  if(reaction_mobile)
262  reaction_mobile->zero_time_step();
263 
265  reaction_immobile->zero_time_step();
266 
267 }
268 
270 {
271  //setting initial condition for immobile concentration matrix
272  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
273  { // Optimize: SWAP LOOPS
274  unsigned int index = el_4_loc_[loc_el];
275  ElementAccessor<3> ele_acc = mesh_->element_accessor(index);
276 
277  for (unsigned int sbi=0; sbi < substances_.size(); sbi++)
278  {
279  conc_immobile[sbi][loc_el] = data_.init_conc_immobile[sbi].value(ele_acc.centre(), ele_acc);
280  }
281  }
282 }
283 
285 {
287 
288  START_TIMER("dual_por_exchange_step");
289  for (unsigned int loc_el = 0; loc_el < distribution_->lsize(); loc_el++)
290  {
292  }
293  END_TIMER("dual_por_exchange_step");
294 
295  if(reaction_mobile) reaction_mobile->update_solution();
296  if(reaction_immobile) reaction_immobile->update_solution();
297 }
298 
299 
300 double **DualPorosity::compute_reaction(double **concentrations, int loc_el)
301 {
302  unsigned int sbi;
303  double conc_average, // weighted (by porosity) average of concentration
304  conc_mob, conc_immob, // new mobile and immobile concentration
305  previous_conc_mob, previous_conc_immob, // mobile and immobile concentration in previous time step
306  conc_max, //difference between concentration and average concentration
307  por_mob, por_immob; // mobile and immobile porosity
308 
309  // get data from fields
311  por_mob = data_.porosity.value(ele.centre(),ele);
312  por_immob = data_.porosity_immobile.value(ele.centre(),ele);
313  arma::Col<double> diff_vec(substances_.size());
314  for (sbi=0; sbi<substances_.size(); sbi++) // Optimize: SWAP LOOPS
315  diff_vec[sbi] = data_.diffusion_rate_immobile[sbi].value(ele.centre(), ele);
316 
317  // if porosity_immobile == 0 then mobile concentration stays the same
318  // and immobile concentration cannot change
319  if (por_immob == 0.0) return conc_immobile;
320 
321  double exponent,
322  temp_exponent = (por_mob + por_immob) / (por_mob * por_immob) * time_->dt();
323 
324  for (sbi = 0; sbi < substances_.size(); sbi++) //over all substances
325  {
326  exponent = diff_vec[sbi] * temp_exponent;
327  //previous values
328  previous_conc_mob = concentration_matrix_[sbi][loc_el];
329  previous_conc_immob = conc_immobile[sbi][loc_el];
330 
331  // ---compute average concentration------------------------------------------
332  conc_average = ((por_mob * previous_conc_mob) + (por_immob * previous_conc_immob))
333  / (por_mob + por_immob);
334 
335  conc_max = std::max(previous_conc_mob-conc_average, previous_conc_immob-conc_average);
336 
337  // the following 2 conditions guarantee:
338  // 1) stability of forward Euler's method
339  // 2) the error of forward Euler's method will not be large
340  if(time_->dt() <= por_mob*por_immob/(max(diff_vec)*(por_mob+por_immob)) &&
341  conc_max <= (2*scheme_tolerance_/(exponent*exponent)*conc_average)) // forward euler
342  {
343  double temp = diff_vec[sbi]*(previous_conc_immob - previous_conc_mob) * time_->dt();
344  // ---compute concentration in mobile area
345  conc_mob = temp / por_mob + previous_conc_mob;
346 
347  // ---compute concentration in immobile area
348  conc_immob = -temp / por_immob + previous_conc_immob;
349  }
350  else //analytic solution
351  {
352  double temp = exp(-exponent);
353  // ---compute concentration in mobile area
354  conc_mob = (previous_conc_mob - conc_average) * temp + conc_average;
355 
356  // ---compute concentration in immobile area
357  conc_immob = (previous_conc_immob - conc_average) * temp + conc_average;
358  }
359 
360  concentration_matrix_[sbi][loc_el] = conc_mob;
361  conc_immobile[sbi][loc_el] = conc_immob;
362  }
363 
364  return conc_immobile;
365 }
366 
367 
369 {
370  int sbi, n_subst;
371  n_subst = substances_.size();
372 
373  vconc_immobile = new Vec [n_subst];
374 
375 
376  for (sbi = 0; sbi < n_subst; sbi++) {
377  VecCreateMPIWithArray(PETSC_COMM_WORLD,1, distribution_->lsize(), mesh_->n_elements(), conc_immobile[sbi],
378  &vconc_immobile[sbi]);
379  VecZeroEntries(vconc_immobile[sbi]);
380 
381  // if(rank == 0)
382  VecZeroEntries(conc_immobile_out[sbi].petsc_vec());
383  }
384 
385  // create output vector scatter
386  IS is;
387  ISCreateGeneral(PETSC_COMM_SELF, mesh_->n_elements(), row_4_el_, PETSC_COPY_VALUES, &is); //WithArray
388  VecScatterCreate(vconc_immobile[0], is, conc_immobile_out[0].petsc_vec(), PETSC_NULL, &vconc_out_scatter);
389  ISDestroy(&(is));
390 }
391 
392 
394 {
395  unsigned int sbi;
396 
397  for (sbi = 0; sbi < substances_.size(); sbi++) {
398  VecScatterBegin(vconc_out_scatter, vconc_immobile[sbi], conc_immobile_out[sbi].petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
399  VecScatterEnd(vconc_out_scatter, vconc_immobile[sbi], conc_immobile_out[sbi].petsc_vec(), INSERT_VALUES, SCATTER_FORWARD);
400  }
401 }
402 
403 
405 {
409  }
410 
411  for (unsigned int sbi = 0; sbi < substances_.size(); sbi++) {
413  }
414 
415  // Register fresh output data
417 
418  if (time_->tlevel() !=0) {
419  // zero_time_step call zero_time_Step of subreactions which performs its own output
420  if (reaction_mobile) reaction_mobile->output_data();
421  if (reaction_immobile) reaction_immobile->output_data();
422  }
423 }
424 
425 
426 bool DualPorosity::evaluate_time_constraint(double &time_constraint)
427 {
428  bool cfl_changed = false;
429  if (reaction_mobile)
430  {
431  if (reaction_mobile->evaluate_time_constraint(time_constraint))
432  cfl_changed = true;
433  }
434  if (reaction_immobile)
435  {
436  double cfl_immobile;
437  if (reaction_immobile->evaluate_time_constraint(cfl_immobile))
438  {
439  time_constraint = std::min(time_constraint, cfl_immobile);
440  cfl_changed = true;
441  }
442  }
443 
444  return cfl_changed;
445 }
TimeGovernor & time()
Definition: equation.hh:148
void output_type(OutputTime::DiscreteSpace rt)
Definition: field_set.hh:211
unsigned int size() const
get global size
Sorption model in immobile zone in case dual porosity is considered.
Definition: sorption.hh:132
FieldSet * eq_data_
Definition: equation.hh:232
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.
std::vector< VectorMPI > conc_immobile_out
concentration array output for immobile phase (gathered - sequential)
Accessor to input data conforming to declared Array.
Definition: accessors.hh:567
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:598
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:61
void output(TimeStep step)
Abstract linear system class.
Definition: balance.hh:35
SubstanceList substances_
Names belonging to substances.
EquationOutput output_fields
Fields indended for output, i.e. all input fields plus those representing solution.
void zero_time_step() override
Definition: mesh.h:80
Field< 3, FieldValue< 3 >::Scalar > porosity_immobile
Immobile porosity field.
Iterator< Ret > find(const string &key) const
~DualPorosity(void)
Destructor.
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
void allocate_output_mpi(void)
Allocates petsc vectors, prepares them for output and creates output vector scatter.
const TimeStep & step(int index=-1) const
static Input::Type::Abstract & it_abstract_mobile_term()
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: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.
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.
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
Definition: accessors.hh:285
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.
void initialize() override
Prepares the object to usage.
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
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:292
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.
std::vector< std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > > output_field_ptr
Fields correspond with conc_immobile_out.
void initialize_fields()
Initializes field sets.
static Input::Type::Abstract & it_abstract_immobile_term()
Mesh * mesh_
Definition: equation.hh:223
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:389
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
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:459
double ** conc_immobile
static Input::Type::Abstract & it_abstract_term()
Support classes for parallel programing.
Sorption model in mobile zone in case dual porosity is considered.
Definition: sorption.hh:106
void output_vector_gather(void) override
Gathers all the parallel vectors to enable them to be output.
FieldCommon & description(const string &description)
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
void initialize(std::shared_ptr< OutputTime > stream, Mesh *mesh, Input::Record in_rec, const TimeGovernor &tg)
Input::Record input_record_
Definition: equation.hh:225
void set_components(const std::vector< string > &names)
Definition: field_set.hh:177
double dt() const
LongIdx * row_4_el_
Indices of rows belonging to elements.
double ** compute_reaction(double **concentrations, int loc_el) override
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:157
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
FieldCommon & name(const string &name)
#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.
static const int registrar
Registrar of class to factory.
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:183
MultiField< 3, FieldValue< 3 >::Scalar > init_conc_immobile
Initial concentrations in the immobile zone.
Record type proxy class.
Definition: type_record.hh:182
bool is_field_output_time(const FieldCommon &field, TimeStep step) const
FieldCommon & set_limits(double min, double max=std::numeric_limits< double >::max())
Field< 3, FieldValue< 3 >::Scalar > porosity
Porosity field.
FieldCommon & flags(FieldFlag::Flags::Mask mask)
LongIdx * el_4_loc_
Indices of elements belonging to local dofs.
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:55
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:96
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:224
FieldSet input_data_set_
unsigned int lsize(int proc) const
get local size