Flow123d  release_3.0.0-1027-g6cabdfa
hc_explicit_sequential.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 hc_explicit_sequential.cc
15  * @brief
16  * @author Jan Brezina
17  */
18 
21 //#include "flow/darcy_flow_mh_output.hh"
22 // TODO:
23 // After having general default values:
24 // make TransportNoting default for AdvectionProcessBase abstract
25 // use default "{}" for secondary equation.
26 // Then we can remove following include.
27 
29 #include "fields/field_common.hh"
30 #include "transport/heat_model.hh"
31 
32 #include "fields/field_set.hh"
33 #include "mesh/mesh.h"
34 #include "io/msh_gmshreader.h"
35 #include "system/sys_profiler.hh"
36 #include "input/input_type.hh"
37 #include "input/accessors.hh"
38 
39 
40 FLOW123D_FORCE_LINK_IN_PARENT(transportOperatorSplitting);
41 FLOW123D_FORCE_LINK_IN_PARENT(concentrationTransportModel);
42 FLOW123D_FORCE_LINK_IN_PARENT(convectionTransport);
44 
45 FLOW123D_FORCE_LINK_IN_PARENT(darcy_flow_mh);
46 FLOW123D_FORCE_LINK_IN_PARENT(richards_lmh);
47 FLOW123D_FORCE_LINK_IN_PARENT(coupling_iterative);
48 
49 
50 namespace it = Input::Type;
51 
53  return it::Abstract("Coupling_Base", "The root record of description of particular the problem to solve.")
54  .close();
55 }
56 
57 
59  return it::Record("Coupling_Sequential",
60  "Record with data for a general sequential coupling.\n")
62  .declare_key("description",it::String(),
63  "Short description of the solved problem.\n"
64  "Is displayed in the main log, and possibly in other text output files.")
66  "Computational mesh common to all equations.")
68  "Simulation time frame and time step.")
71  "Flow equation, provides the velocity field as a result.")
73  "Transport of soluted substances, depends on the velocity field from a Flow equation.")
75  "Heat transfer, depends on the velocity field from a Flow equation.")
76  .close();
77 }
78 
79 
81 
82 
83 
84 std::shared_ptr<AdvectionProcessBase> HC_ExplicitSequential::make_advection_process(string process_key)
85 {
86  using namespace Input;
87  // setup heat object
88  Iterator<AbstractRecord> it = in_record_.find<AbstractRecord>(process_key);
89 
90  if (it) {
91  auto process = (*it).factory< AdvectionProcessBase, Mesh &, const Input::Record >(*mesh, *it);
92 
93  // setup fields
94  process->data()["cross_section"]
95  .copy_from(water->data()["cross_section"]);
96  /*
97  if (water_content_saturated_) // only for unsteady Richards water model
98  process->data()["porosity"].copy_from(*water_content_saturated_);
99 
100  if (water_content_p0_)
101  process->data()["water_content"].copy_from(*water_content_p0_);
102  else {
103 
104  }*/
105 
106  FieldCommon *porosity = process->data().field("porosity");
107  process->data()["water_content"].copy_from( *porosity );
108 
109 
110  process->initialize();
111  return process;
112  } else {
113  return std::make_shared<TransportNothing>(*mesh);
114  }
115 }
116 
117 /**
118  * FUNCTION "MAIN" FOR COMPUTING MIXED-HYBRID PROBLEM FOR UNSTEADY SATURATED FLOW
119  */
121 : in_record_(in_record)
122 {
123  START_TIMER("HC constructor");
124  using namespace Input;
125 
126  // Read mesh
127  {
128  START_TIMER("HC read mesh");
129 
130  mesh = BaseMeshReader::mesh_factory( in_record.val<Record>("mesh") );
131 
132  //getting description for the Profiler
133  string description;
134  in_record.opt_val<string>("description", description);
135 
137  description,
138  //"Description has to be set in main. by different method.",
139  mesh->n_elements());
140  }
141 
142  // setup primary equation - water flow object
143  AbstractRecord prim_eq = in_record.val<AbstractRecord>("flow_equation");
144  // Need explicit template types here, since reference is used (automatically passing by value)
145  water = prim_eq.factory< DarcyFlowInterface, Mesh &, const Input::Record>(*mesh, prim_eq);
146  water->initialize();
147  std::stringstream ss; // print warning message with table of uninitialized fields
148  if ( FieldCommon::print_message_table(ss, "HC explicit sequential") ) {
149  WarningOut() << ss.str();
150  }
151 
152  RegionSet bulk_set = mesh->region_db().get_region_set("BULK");
153  water_content_saturated_ = water->data().field("water_content_saturated");
155  water_content_saturated_ = nullptr;
156 
157  water_content_p0_ = water->data().field("water_content_p0");
159  water_content_p0_ = nullptr;
160 
161  processes_.push_back(AdvectionData(make_advection_process("solute_equation")));
162  processes_.push_back(AdvectionData(make_advection_process("heat_equation")));
163 }
164 
166 {
167  if (pdata.process->time().is_end()) return;
168 
169  is_end_all_=false;
170  if ( pdata.process->time().step().le(pdata.velocity_time) ) {
171  // having information about velocity field we can perform transport step
172 
173  // here should be interpolation of the velocity at least if the interpolation time
174  // is not close to the solved_time of the water module
175  // for simplicity we use only last velocity field
176  if (pdata.velocity_changed) {
177  //DBGMSG("velocity update\n");
178  pdata.process->set_velocity_field( water->get_mh_dofhandler() );
179  pdata.velocity_changed = false;
180  }
181  if (pdata.process->time().tlevel() == 0) pdata.process->zero_time_step();
182 
183  pdata.process->update_solution();
184  pdata.process->output_data();
185  }
186 }
187 
188 
189 /**
190  * TODO:
191  * - have support for steady problems in TimeGovernor, make Noting problems steady
192  * - apply splitting of compute_one_step to particular models
193  * - how to set output time marks for steady problems (we need solved time == infinity) but
194  * add no time marks
195  * - allow create steady time governor without time marks (at least in nothing models)
196  * - pass refference to time marks in time governor constructor?
197  */
198 
200 {
201  START_TIMER("HC run simulation");
202  // following should be specified in constructor:
203  // value for velocity interpolation :
204  // theta = 0 velocity from beginning of transport interval (fully explicit method)
205  // theta = 0.5 velocity from center of transport interval ( mimic Crank-Nicholson)
206  // theta = 1.0 velocity from end of transport interval (partialy explicit scheme)
207  const double theta=0.5;
208 
209  {
210  START_TIMER("HC water zero time step");
211  water->zero_time_step();
212  }
213 
214 
215  // following cycle is designed to support independent time stepping of
216  // both processes. The question is which value of the water field use to compute a transport step.
217  // Meaningful cases are
218  // 1) beginning (fully explicit method)
219  // 2) center ( mimic Crank-Nicholson)
220  // 3) end of the interval (partialy explicit scheme)
221  // However with current implementation of the explicit transport on have to assembly transport matrix for
222  // every new value of the velocity field. So we have to keep same velocity field over some time interval t_dt
223  // which is further split into shorter time intervals ts_dt dictated by the CFL condition.
224  // One can consider t_dt as the transport time step and apply one of the previous three cases.
225  //
226  // The question is how to choose intervals t_dt. That should depend on variability of the velocity field in time.
227  // Currently we simply use t_dt == w_dt.
228 
229  is_end_all_=false;
230  while (! is_end_all_) {
231  is_end_all_ = true;
232 
233  double water_dt=water->time().estimate_dt();
234  if (water->time().is_end()) water_dt = TimeGovernor::inf_time;
235 
236  // in future here could be re-estimation of transport planed time according to
237  // evolution of the velocity field. Consider the case w_dt << t_dt and velocity almost constant in time
238  // which suddenly rise in time 3*w_dt. First we the planed transport time step t_dt could be quite big, but
239  // in time 3*w_dt we can reconsider value of t_dt to better capture changing velocity.
241  for(auto &pdata : processes_) {
242  pdata.process->set_time_upper_constraint(water_dt, "Flow time step");
243  pdata.velocity_time = theta * pdata.process->planned_time() + (1-theta) * pdata.process->solved_time();
244  min_velocity_time = min(min_velocity_time, pdata.velocity_time);
245  }
246 
247  // printing water and transport times every step
248  //MessageOut().fmt("HC_EXPL_SEQ: velocity_interpolation_time: {}, water_time: {} transport time: {}\n",
249  // velocity_interpolation_time, water->time().t(), transport_reaction->time().t());
250 
251  // if transport is off, transport should return infinity solved and planned times so that
252  // only water branch takes the place
253 
254  if (! water->time().is_end() ) {
255  is_end_all_=false;
256  if (water->solved_time() < min_velocity_time) {
257 
258  // solve water over the nearest transport interval
259  water->update_solution();
260 
261  // here possibly save solution from water for interpolation in time
262 
263  //water->time().view("WATER"); //show water time governor
264 
265  //water->output_data();
266  water->choose_next_time();
267 
268  for(auto &process : processes_) process.velocity_changed = true;
269 
270  }
271  }
272  advection_process_step(processes_[0]); // solute
273  advection_process_step(processes_[1]); // heat
274  }
275  //MessageOut().fmt("End of simulation at time: {}\n", max(solute->solved_time(), heat->solved_time()));
276 }
277 
278 
280  water.reset();
281  for(auto &pdata : processes_) pdata.process.reset();
282  delete mesh;
283 }
284 
285 
286 
287 
288 //-----------------------------------------------------------------------------
289 // vim: set cindent:
290 //-----------------------------------------------------------------------------
Mesh * mesh
mesh common to darcy flow and transport
Common abstract parent of all Field<...> classes.
Definition: field_common.hh:73
virtual void copy_from(const FieldCommon &other)=0
static Input::Type::Abstract & get_input_type()
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:598
Abstract linear system class.
Definition: balance.hh:35
static const Input::Type::Record & get_input_type()
Definition: mesh.cc:76
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
Definition: mesh.h:76
Abstract & close()
Close the Abstract and add its to type repository (see TypeRepository::add_type). ...
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:303
HC_ExplicitSequential(Input::Record in_record)
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:195
static Default optional()
The factory function to make an empty default value which is optional.
Definition: type_record.hh:124
bool opt_val(const string &key, Ret &value) const
static Input::Type::Abstract & get_input_type()
Common specification of the input record for secondary equations.
static Mesh * mesh_factory(const Input::Record &input_mesh_rec)
AdvectionPtr make_advection_process(std::string process_key)
std::vector< AdvectionData > processes_
solute transport with chemistry through operator splitting
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
std::shared_ptr< DarcyFlowInterface > water
steady or unsteady water flow simulator based on MH scheme
const Ret val(const string &key) const
#define START_TIMER(tag)
Starts a timer with specified tag.
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
static Input::Type::Abstract & get_input_type()
FieldCommon * water_content_saturated_
Class for declaration of polymorphic Record.
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:459
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
Definition: mesh.h:346
virtual FieldResult field_result(RegionSet region_set) const =0
Indicates special field states.
static Profiler * instance()
static const Input::Type::Record & get_input_type()
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:246
Discontinuous Galerkin method for equation of transport with dispersion.
static const Input::Type::Record & get_input_type()
void advection_process_step(AdvectionData &pdata)
Record type proxy class.
Definition: type_record.hh:182
const std::shared_ptr< Type > factory(Arguments...arguments) const
FLOW123D_FORCE_LINK_IN_PARENT(transportOperatorSplitting)
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:96
Class for declaration of the input data that are in string format.
Definition: type_base.hh:589
static const double inf_time
Infinity time used for steady case.
void set_task_info(string description, int size)