Flow123d  JS_before_hm-919-g5f1bbbf
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 #include "system/index_types.hh"
24 
26 #include <petscmat.h>
27 
28 #include "io/output_time.hh"
29 #include "tools/time_governor.hh"
30 #include "coupling/equation.hh"
31 #include "coupling/balance.hh"
32 #include "transport/transport.h"
33 #include "mesh/mesh.h"
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\"")
77  .declare_key("balance", Balance::get_input_type(), Default("{}"),
78  "Settings for computing mass balance.")
79  .declare_key("output_stream", OutputTime::get_input_type(), Default("{}"),
80  "Output stream settings.\n Specify file format, precision etc.")
81  .declare_key("substances", Array( Substance::get_input_type(), 1), Default::obligatory(),
82  "Specification of transported substances.")
83  // input data
84  .declare_key("transport", ConcentrationTransportBase::get_input_type(), Default::obligatory(),
85  "Type of the numerical method for the transport equation.")
86  .declare_key("reaction_term", ReactionTerm::it_abstract_term(), Default::optional(),
87  "Reaction model involved in transport.")
88 /*
89  .declare_key("output_fields", Array(ConvectionTransport::get_output_selection()),
90  Default("\"conc\""),
91  "List of fields to write to output file.")*/
92  .close();
93 }
94 
95 
97  Input::register_class< TransportOperatorSplitting, Mesh &, const Input::Record>("Coupling_OperatorSplitting") +
99 
100 
101 
102 
103 
104 
105 
107 {
108  *this += porosity.name("porosity")
109  .description("Porosity of the mobile phase.")
110  .input_default("1.0")
111  .units( UnitSI::dimensionless() )
112  .flags_add(in_main_matrix & in_rhs);
113 
114  *this += water_content.name("water_content")
115  .description("INTERNAL. Water content passed from unsaturated Darcy flow model.")
116  .units( UnitSI::dimensionless() )
117  .flags_add(input_copy & in_time_term & in_main_matrix & in_rhs);
118 
119  *this += cross_section.name("cross_section")
120  .flags( FieldFlag::input_copy )
121  .flags_add(in_time_term & in_main_matrix & in_rhs);
122 
123  *this += sources_density.name("sources_density")
124  .description("Density of concentration sources.")
125  .input_default("0.0")
126  .units( UnitSI().kg().m(-3).s(-1) )
127  .flags_add(in_rhs);
128 
129  *this += sources_sigma.name("sources_sigma")
130  .description("Concentration flux.")
131  .input_default("0.0")
132  .units( UnitSI().s(-1) )
133  .flags_add(in_main_matrix & in_rhs);
134 
135  *this += sources_conc.name("sources_conc")
136  .description("Concentration sources threshold.")
137  .input_default("0.0")
138  .units( UnitSI().kg().m(-3) )
139  .flags_add(in_rhs);
140 }
141 
142 
143 
144 
145 
146 
147 
148 /////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
149 
150 
151 
153 : AdvectionProcessBase(init_mesh, in_rec),
154 // Semchem_reactions(NULL),
155  cfl_convection(numeric_limits<double>::max()),
156  cfl_reaction(numeric_limits<double>::max())
157 {
158  START_TIMER("TransportOperatorSpliting");
159 
160  Distribution *el_distribution;
161  LongIdx *el_4_loc;
162 
163  Input::AbstractRecord trans = in_rec.val<Input::AbstractRecord>("transport");
164  convection = trans.factory< ConcentrationTransportBase, Mesh &, const Input::Record >(init_mesh, trans);
165 
166  time_ = new TimeGovernor(in_rec.val<Input::Record>("time"), TimeMark::none_type, false);
167  ASSERT( time_->is_default() == false ).error("Missing key 'time' in Coupling_OperatorSplitting.");
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 MixedPtr<FE_P_disc> fe(0);
203  shared_ptr<DOFHandlerMultiDim> dof_handler = make_shared<DOFHandlerMultiDim>(*mesh_);
204  shared_ptr<DiscreteSpace> ds = make_shared<EqualOrderDiscreteSpace>( mesh_, fe);
205  dof_handler->distribute_dofs(ds);
206 
207  reaction->substances(convection->substances())
208  .concentration_matrix(convection->get_concentration_matrix(),
209  el_distribution, el_4_loc, convection->get_row_4_el())
210  .output_stream(convection->output_stream())
211  .set_dh(dof_handler)
212  .set_time_governor((TimeGovernor &)convection->time());
213 
214  reaction->initialize();
215 
216  } else {
217  reaction = nullptr;
218  //Semchem_reactions = nullptr;
219  }
220 }
221 
223 {
224  //delete field_output;
225  //if (Semchem_reactions) delete Semchem_reactions;
226  delete time_;
227 }
228 
229 
231 {
232  //coupling - passing fields
233  if(reaction)
234  if( typeid(*reaction) == typeid(SorptionSimple) ||
235  typeid(*reaction) == typeid(DualPorosity)
236  )
237  {
238  reaction->data().set_field("porosity", convection->data()["porosity"]);
239  }
240 }
241 
242 
243 
244 
246 
247 
248  START_TIMER("TOS-output data");
249 
250 
251  convection->output_data();
252  if(reaction) reaction->output_data(); // do not perform write_time_frame
253  convection->output_stream()->write_time_frame();
254 
255  END_TIMER("TOS-output data");
256 }
257 
258 
260 {
261  //DebugOut() << "tos ZERO TIME STEP.\n";
262  convection->zero_time_step();
263  convection->calculate_concentration_matrix(); // due to reading of init_conc in reactions
264  if(reaction)
265  {
266  reaction->zero_time_step();
267  reaction->output_data(); // do not perform write_time_frame
268  }
269  convection->output_stream()->write_time_frame();
270 
271 }
272 
273 
274 
276 
277  vector<double> source(convection->n_substances()), region_mass(mesh_->region_db().bulk_size());
278 
279  time_->next_time();
280  time_->view("TOS"); //show time governor
281 
282  convection->set_target_time(time_->t());
283 
284  START_TIMER("TOS-one step");
285  int steps=0;
286  while ( convection->time().step().lt(time_->t()) )
287  {
288  steps++;
289  // one internal step
290  // we call evaluate_time_constraint() of convection and reaction separately to
291  // make sure that both routines are executed.
292  bool cfl_convection_changed = convection->evaluate_time_constraint(cfl_convection);
293  bool cfl_reaction_changed = (reaction?reaction->evaluate_time_constraint(cfl_reaction):0);
294  bool cfl_changed = cfl_convection_changed || cfl_reaction_changed;
295 
296  if (steps == 1 || cfl_changed)
297  {
298  convection->time().set_upper_constraint(cfl_convection, "Time step constrained by transport CFL condition (including both flow and sources).");
299  convection->time().set_upper_constraint(cfl_reaction, "Time step constrained by reaction CFL condition.");
300 
301  // fix step with new constraint
302  convection->time().fix_dt_until_mark();
303 
304  convection->time().view("Convection"); // write TG only once on change
305  }
306 
307  convection->update_solution();
308 
309 
310  if (balance_->cumulative())
311  {
312  START_TIMER("TOS-balance");
313 
314  // save mass after transport step
315  for (unsigned int sbi=0; sbi<convection->n_substances(); sbi++)
316  {
317  balance_->calculate_mass(convection->get_subst_idx()[sbi], convection->get_solution(sbi), region_mass);
318  source[sbi] = 0;
319  for (unsigned int ri=0; ri<mesh_->region_db().bulk_size(); ri++)
320  source[sbi] -= region_mass[ri];
321  }
322 
323  END_TIMER("TOS-balance");
324  }
325 
326  if(reaction) {
327  convection->calculate_concentration_matrix();
328  reaction->update_solution();
329  convection->update_after_reactions(true);
330  }
331  else
332  convection->update_after_reactions(false);
333 
334  //if(Semchem_reactions) Semchem_reactions->update_solution();
335 
336 
337 
338  if (balance_->cumulative())
339  {
340  START_TIMER("TOS-balance");
341 
342  for (unsigned int sbi=0; sbi<convection->n_substances(); sbi++)
343  {
344  // compute mass difference due to reactions
345  balance_->calculate_mass(convection->get_subst_idx()[sbi], convection->get_solution(sbi), region_mass);
346  for (unsigned int ri=0; ri<mesh_->region_db().bulk_size(); ri++)
347  source[sbi] += region_mass[ri];
348 
349  // update balance of sources due to reactions
350  balance_->add_cumulative_source(sbi, source[sbi]);
351  }
352 
353  END_TIMER("TOS-balance");
354  }
355  }
356 
357  LogOut().fmt("CONVECTION: steps: {}\n", steps);
358 }
359 
360 
361 
362 
363 
365 {
366  convection->set_velocity_field( flux_field );
367 }
368 
369 
370 
371 
TimeGovernor & time()
Definition: equation.hh:151
FieldSet * eq_data_
Definition: equation.hh:229
static const Input::Type::Record & get_input_type()
Main balance input record type.
Definition: balance.cc:50
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:602
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:40
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:78
Iterator< Ret > find(const string &key) const
static string subfields_address()
const RegionDB & region_db() const
Definition: mesh.h:143
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
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:261
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.
static Input::Type::Record & record_template()
Template Record with common keys for derived equations.
Definition: equation.cc:36
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
virtual void set_velocity_field(std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed >> flux_field) override
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:220
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:232
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
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
static Input::Type::Abstract & it_abstract_term()
Support classes for parallel programing.
void output_data() override
Write computed fields.
std::shared_ptr< ConcentrationTransportBase > convection
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
#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:221
Definition: field.hh:60
virtual ~TransportOperatorSplitting()
Destructor.
static const int registrar
Registrar of class to factory.