Flow123d  release_3.0.0-1264-g45bfb2a
transport_operator_splitting.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 transport_operator_splitting.cc
15  * @brief
16  */
17 
18 #include <iostream>
19 #include <iomanip>
20 
21 #include "system/system.hh"
22 #include "system/sys_profiler.hh"
23 
25 #include <petscmat.h>
26 
27 #include "io/output_time.hh"
28 #include "tools/time_governor.hh"
29 #include "coupling/equation.hh"
30 #include "coupling/balance.hh"
31 #include "transport/transport.h"
32 #include "mesh/mesh.h"
33 #include "mesh/long_idx.hh"
34 
38 #include "reaction/sorption.hh"
40 
41 #include "la/distribution.hh"
42 #include "input/input_type.hh"
43 #include "input/accessors.hh"
44 #include "input/factory.hh"
46 #include "fem/fe_p.hh"
47 
48 FLOW123D_FORCE_LINK_IN_CHILD(transportOperatorSplitting)
49 
50 FLOW123D_FORCE_LINK_IN_PARENT(firstOrderReaction)
51 FLOW123D_FORCE_LINK_IN_PARENT(radioactiveDecay)
53 FLOW123D_FORCE_LINK_IN_PARENT(sorptionMobile)
54 FLOW123D_FORCE_LINK_IN_PARENT(sorptionImmobile)
56 
57 
58 using namespace Input::Type;
59 
60 
61 
62 Abstract & ConcentrationTransportBase::get_input_type() {
63  return Abstract("Solute",
64  "Transport of soluted substances.")
65  .close();
66 }
67 
68 
70  return Record("Coupling_OperatorSplitting",
71  "Transport by convection and/or diffusion\n"
72  "coupled with reaction and adsorption model (ODE per element)\n"
73  " via operator splitting.")
75  .add_attribute( FlowAttribute::subfields_address(), "\"/problem/solute_equation/substances/*/name\"")
76  .declare_key("time", TimeGovernor::get_input_type(), Default::obligatory(),
77  "Time governor settings for the transport equation.")
78  .declare_key("balance", Balance::get_input_type(), Default("{}"),
79  "Settings for computing mass balance.")
80  .declare_key("output_stream", OutputTime::get_input_type(), Default("{}"),
81  "Output stream settings.\n Specify file format, precision etc.")
82  .declare_key("substances", Array( Substance::get_input_type(), 1), Default::obligatory(),
83  "Specification of transported substances.")
84  // input data
85  .declare_key("transport", ConcentrationTransportBase::get_input_type(), Default::obligatory(),
86  "Type of the numerical method for the transport equation.")
87  .declare_key("reaction_term", ReactionTerm::it_abstract_term(), Default::optional(),
88  "Reaction model involved in transport.")
89 /*
90  .declare_key("output_fields", Array(ConvectionTransport::get_output_selection()),
91  Default("\"conc\""),
92  "List of fields to write to output file.")*/
93  .close();
94 }
95 
96 
98  Input::register_class< TransportOperatorSplitting, Mesh &, const Input::Record>("Coupling_OperatorSplitting") +
100 
101 
102 
103 
104 
105 
106 
108 {
109  *this += porosity.name("porosity")
110  .description("Porosity of the mobile phase.")
111  .input_default("1.0")
112  .units( UnitSI::dimensionless() )
113  .flags_add(in_main_matrix & in_rhs);
114 
115  *this += water_content.name("water_content")
116  .description("INTERNAL. Water content passed from unsaturated Darcy flow model.")
117  .units( UnitSI::dimensionless() )
118  .flags_add(input_copy & in_time_term & in_main_matrix & in_rhs);
119 
120  *this += cross_section.name("cross_section")
121  .flags( FieldFlag::input_copy )
122  .flags_add(in_time_term & in_main_matrix & in_rhs);
123 
124  *this += sources_density.name("sources_density")
125  .description("Density of concentration sources.")
126  .input_default("0.0")
127  .units( UnitSI().kg().m(-3).s(-1) )
128  .flags_add(in_rhs);
129 
130  *this += sources_sigma.name("sources_sigma")
131  .description("Concentration flux.")
132  .input_default("0.0")
133  .units( UnitSI().s(-1) )
134  .flags_add(in_main_matrix & in_rhs);
135 
136  *this += sources_conc.name("sources_conc")
137  .description("Concentration sources threshold.")
138  .input_default("0.0")
139  .units( UnitSI().kg().m(-3) )
140  .flags_add(in_rhs);
141 }
142 
143 
144 
145 
146 
147 
148 
149 /////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
150 
151 
152 
154 : AdvectionProcessBase(init_mesh, in_rec),
155 // Semchem_reactions(NULL),
156  cfl_convection(numeric_limits<double>::max()),
157  cfl_reaction(numeric_limits<double>::max())
158 {
159  START_TIMER("TransportOperatorSpliting");
160 
161  Distribution *el_distribution;
162  LongIdx *el_4_loc;
163 
164  Input::AbstractRecord trans = in_rec.val<Input::AbstractRecord>("transport");
165  convection = trans.factory< ConcentrationTransportBase, Mesh &, const Input::Record >(init_mesh, trans);
166 
167  time_ = new TimeGovernor(in_rec.val<Input::Record>("time"), TimeMark::none_type, false);
168  convection->set_time_governor(time());
169 
170  // Initialize list of substances.
171  convection->substances().initialize(in_rec.val<Input::Array>("substances"));
172 
173  // Initialize output stream.
174  convection->set_output_stream(OutputTime::create_output_stream("solute", in_rec.val<Input::Record>("output_stream"), time().get_unit_string()));
175 
176 
177  // initialization of balance object
178 
179  balance_ = make_shared<Balance>("mass", mesh_);
180  balance_->init_from_input(in_rec.val<Input::Record>("balance"), this->time());
181 
182  balance_->units(UnitSI().kg(1));
183  convection->set_balance_object(balance_);
184 
185  convection->initialize(); //
186 
187 
188 
189 
190  time_ = new TimeGovernor(in_rec.val<Input::Record>("time"), convection->mark_type());
191 
192  this->eq_data_ = &(convection->data());
193 
194  convection->get_par_info(el_4_loc, el_distribution);
195  Input::Iterator<Input::AbstractRecord> reactions_it = in_rec.find<Input::AbstractRecord>("reaction_term");
196  if ( reactions_it ) {
197  // TODO: allowed instances in this case are only
198  // FirstOrderReaction, RadioactiveDecay, SorptionSimple and DualPorosity
199  reaction = (*reactions_it).factory< ReactionTerm, Mesh &, Input::Record >(init_mesh, *reactions_it);
200 
201  //initialization of DOF handler
202  static FE_P_disc<0> fe0(0);
203  static FE_P_disc<1> fe1(0);
204  static FE_P_disc<2> fe2(0);
205  static FE_P_disc<3> fe3(0);
206  shared_ptr<DOFHandlerMultiDim> dof_handler = make_shared<DOFHandlerMultiDim>(*mesh_);
207  shared_ptr<DiscreteSpace> ds = make_shared<EqualOrderDiscreteSpace>( mesh_, &fe0, &fe1, &fe2, &fe3);
208  dof_handler->distribute_dofs(ds);
209 
210  reaction->substances(convection->substances())
211  .concentration_matrix(convection->get_concentration_matrix(),
212  el_distribution, el_4_loc, convection->get_row_4_el())
213  .output_stream(convection->output_stream())
214  .set_dh(dof_handler)
215  .set_time_governor((TimeGovernor &)convection->time());
216 
217  reaction->initialize();
218 
219  } else {
220  reaction = nullptr;
221  //Semchem_reactions = nullptr;
222  }
223 }
224 
226 {
227  //delete field_output;
228  //if (Semchem_reactions) delete Semchem_reactions;
229  delete time_;
230 }
231 
232 
234 {
235  //coupling - passing fields
236  if(reaction)
237  if( typeid(*reaction) == typeid(SorptionSimple) ||
238  typeid(*reaction) == typeid(DualPorosity)
239  )
240  {
241  reaction->data().set_field("porosity", convection->data()["porosity"]);
242  }
243 }
244 
245 
246 
247 
249 
250 
251  START_TIMER("TOS-output data");
252 
253 
254  convection->output_data();
255  if(reaction) reaction->output_data(); // do not perform write_time_frame
256  convection->output_stream()->write_time_frame();
257 
258  END_TIMER("TOS-output data");
259 }
260 
261 
263 {
264  //DebugOut() << "tos ZERO TIME STEP.\n";
265  convection->zero_time_step();
266  convection->calculate_concentration_matrix(); // due to reading of init_conc in reactions
267  if(reaction)
268  {
269  reaction->zero_time_step();
270  reaction->output_data(); // do not perform write_time_frame
271  }
272  convection->output_stream()->write_time_frame();
273 
274 }
275 
276 
277 
279 
280  vector<double> source(convection->n_substances()), region_mass(mesh_->region_db().bulk_size());
281 
282  time_->next_time();
283  time_->view("TOS"); //show time governor
284 
285  convection->set_target_time(time_->t());
286 
287  START_TIMER("TOS-one step");
288  int steps=0;
289  while ( convection->time().step().lt(time_->t()) )
290  {
291  steps++;
292  // one internal step
293  // we call evaluate_time_constraint() of convection and reaction separately to
294  // make sure that both routines are executed.
295  bool cfl_convection_changed = convection->evaluate_time_constraint(cfl_convection);
296  bool cfl_reaction_changed = (reaction?reaction->evaluate_time_constraint(cfl_reaction):0);
297  bool cfl_changed = cfl_convection_changed || cfl_reaction_changed;
298 
299  if (steps == 1 || cfl_changed)
300  {
301  convection->time().set_upper_constraint(cfl_convection, "Time step constrained by transport CFL condition (including both flow and sources).");
302  convection->time().set_upper_constraint(cfl_reaction, "Time step constrained by reaction CFL condition.");
303 
304  // fix step with new constraint
305  convection->time().fix_dt_until_mark();
306 
307  convection->time().view("Convection"); // write TG only once on change
308  }
309 
310  convection->update_solution();
311 
312 
313  if (balance_->cumulative())
314  {
315  START_TIMER("TOS-balance");
316 
317  // save mass after transport step
318  for (unsigned int sbi=0; sbi<convection->n_substances(); sbi++)
319  {
320  balance_->calculate_mass(convection->get_subst_idx()[sbi], convection->get_solution(sbi), region_mass);
321  source[sbi] = 0;
322  for (unsigned int ri=0; ri<mesh_->region_db().bulk_size(); ri++)
323  source[sbi] -= region_mass[ri];
324  }
325 
326  END_TIMER("TOS-balance");
327  }
328 
329  if(reaction) {
330  convection->calculate_concentration_matrix();
331  reaction->update_solution();
332  convection->update_after_reactions(true);
333  }
334  else
335  convection->update_after_reactions(false);
336 
337  //if(Semchem_reactions) Semchem_reactions->update_solution();
338 
339 
340 
341  if (balance_->cumulative())
342  {
343  START_TIMER("TOS-balance");
344 
345  for (unsigned int sbi=0; sbi<convection->n_substances(); sbi++)
346  {
347  // compute mass difference due to reactions
348  balance_->calculate_mass(convection->get_subst_idx()[sbi], convection->get_solution(sbi), region_mass);
349  for (unsigned int ri=0; ri<mesh_->region_db().bulk_size(); ri++)
350  source[sbi] += region_mass[ri];
351 
352  // update balance of sources due to reactions
353  balance_->add_cumulative_source(sbi, source[sbi]);
354  }
355 
356  END_TIMER("TOS-balance");
357  }
358  }
359 
360  LogOut().fmt("CONVECTION: steps: {}\n", steps);
361 }
362 
363 
364 
365 
366 
368 {
369  convection->set_velocity_field( dh );
370 }
371 
372 
373 
374 
TimeGovernor & time()
Definition: equation.hh:148
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
FieldSet * eq_data_
Definition: equation.hh:232
static const Input::Type::Record & get_input_type()
Main balance input record type.
Definition: balance.cc:48
Abstract base class for equation clasess.
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
unsigned int bulk_size() const
Definition: region.cc:269
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportOperatorSplittiong.
static const Input::Type::Record & get_input_type()
The specification of output stream.
Definition: output_time.cc:37
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Abstract linear system class.
Definition: balance.hh:37
virtual void set_velocity_field(const MH_DofHandler &dh) override
void next_time()
Proceed to the next time according to current estimated time step.
static std::shared_ptr< OutputTime > create_output_stream(const std::string &equation_name, const Input::Record &in_rec, std::string unit_str)
This method delete all object instances of class OutputTime stored in output_streams vector...
Definition: output_time.cc:184
Definition: mesh.h:76
Iterator< Ret > find(const string &key) const
static string subfields_address()
const RegionDB & region_db() const
Definition: mesh.h:141
static const Input::Type::Record & get_input_type()
Input type for a substance.
Definition: substance.cc:29
TransportOperatorSplitting(Mesh &init_mesh, const Input::Record in_rec)
Constructor.
static Input::Type::Abstract & get_input_type()
Common specification of the input record for secondary equations.
double t() const
Abstract & close()
Close the Abstract and add its to type repository (see TypeRepository::add_type). ...
Basic time management functionality for unsteady (and steady) solvers (class Equation).
static const Type none_type
Mark Type with all bits unset.
Definition: time_marks.hh:95
#define LogOut()
Macro defining &#39;log&#39; record of log.
Definition: logger.hh:249
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
Basic time management class.
Class for declaration of inputs sequences.
Definition: type_base.hh:346
void view(const char *name="") const
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:196
static Input::Type::Abstract & get_input_type()
Common specification of the input record for secondary equations.
std::shared_ptr< ReactionTerm > reaction
Record & add_attribute(std::string key, TypeBase::json_string value)
Add TYPE key as obligatory.
Definition: type_record.cc:354
static constexpr Mask input_copy
Definition: field_flag.hh:44
Definitions of basic Lagrangean finite elements with polynomial shape functions.
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.
Mesh * mesh_
Definition: equation.hh:223
double cfl_convection
Time restriction due to transport.
std::shared_ptr< Balance > balance_
object for calculation and writing the mass balance to file.
Definition: equation.hh:235
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:503
Class for declaration of polymorphic Record.
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
static Input::Type::Abstract & it_abstract_term()
Support classes for parallel programing.
void output_data() override
Write computed fields.
std::shared_ptr< ConcentrationTransportBase > convection
static const Input::Type::Record & get_input_type()
#define END_TIMER(tag)
Ends a timer with specified tag.
double cfl_reaction
Time restriction due to reactions.
Simple sorption model without dual porosity.
Definition: sorption.hh:52
Record type proxy class.
Definition: type_record.hh:182
#define FLOW123D_FORCE_LINK_IN_PARENT(x)
Definition: global_defs.h:183
const std::shared_ptr< Type > factory(Arguments...arguments) const
Class for representation SI units of Fields.
Definition: unit_si.hh:40
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
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
virtual ~TransportOperatorSplitting()
Destructor.
static const int registrar
Registrar of class to factory.