Flow123d  release_2.2.0-33-g759111d
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 
48 
49 namespace it = Input::Type;
50 
52  return it::Abstract("Coupling_Base", "The root record of description of particular the problem to solve.")
53  .close();
54 }
55 
56 
58  return it::Record("Coupling_Sequential",
59  "Record with data for a general sequential coupling.\n")
61  .declare_key("description",it::String(),
62  "Short description of the solved problem.\n"
63  "Is displayed in the main log, and possibly in other text output files.")
65  "Computational mesh common to all equations.")
67  "Simulation time frame and time step.")
70  "Flow equation, provides the velocity field as a result.")
72  "Transport of soluted substances, depends on the velocity field from a Flow equation.")
74  "Heat transfer, depends on the velocity field from a Flow equation.")
75  .close();
76 }
77 
78 
80 
81 
82 
83 std::shared_ptr<AdvectionProcessBase> HC_ExplicitSequential::make_advection_process(string process_key)
84 {
85  using namespace Input;
86  // setup heat object
87  Iterator<AbstractRecord> it = in_record_.find<AbstractRecord>(process_key);
88 
89  if (it) {
90  auto process = (*it).factory< AdvectionProcessBase, Mesh &, const Input::Record >(*mesh, *it);
91 
92  // setup fields
93  process->data()["cross_section"]
94  .copy_from(water->data()["cross_section"]);
95  /*
96  if (water_content_saturated_) // only for unsteady Richards water model
97  process->data()["porosity"].copy_from(*water_content_saturated_);
98 
99  if (water_content_p0_)
100  process->data()["water_content"].copy_from(*water_content_p0_);
101  else {
102 
103  }*/
104 
105  FieldCommon *porosity = process->data().field("porosity");
106  process->data()["water_content"].copy_from( *porosity );
107 
108 
109  process->initialize();
110  return process;
111  } else {
112  return std::make_shared<TransportNothing>(*mesh);
113  }
114 }
115 
116 /**
117  * FUNCTION "MAIN" FOR COMPUTING MIXED-HYBRID PROBLEM FOR UNSTEADY SATURATED FLOW
118  */
120 : in_record_(in_record)
121 {
122  START_TIMER("HC constructor");
123  using namespace Input;
124 
125  // Read mesh
126  {
127  START_TIMER("HC read mesh");
128 
129  mesh = BaseMeshReader::mesh_factory( in_record.val<Record>("mesh") );
130 
131  //getting description for the Profiler
132  string description;
133  in_record.opt_val<string>("description", description);
134 
136  description,
137  //"Description has to be set in main. by different method.",
138  mesh->n_elements());
139  }
140 
141  // setup primary equation - water flow object
142  AbstractRecord prim_eq = in_record.val<AbstractRecord>("flow_equation");
143  // Need explicit template types here, since reference is used (automatically passing by value)
144  water = prim_eq.factory< DarcyFlowInterface, Mesh &, const Input::Record>(*mesh, prim_eq);
145  water->initialize();
146  std::stringstream ss; // print warning message with table of uninitialized fields
147  if ( FieldCommon::print_message_table(ss, "HC explicit sequential") ) {
148  WarningOut() << ss.str();
149  }
150 
151  RegionSet bulk_set = mesh->region_db().get_region_set("BULK");
152  water_content_saturated_ = water->data().field("water_content_saturated");
154  water_content_saturated_ = nullptr;
155 
156  water_content_p0_ = water->data().field("water_content_p0");
158  water_content_p0_ = nullptr;
159 
160  processes_.push_back(AdvectionData(make_advection_process("solute_equation")));
161  processes_.push_back(AdvectionData(make_advection_process("heat_equation")));
162 }
163 
165 {
166  if (pdata.process->time().is_end()) return;
167 
168  is_end_all_=false;
169  if ( pdata.process->time().step().le(pdata.velocity_time) ) {
170  // having information about velocity field we can perform transport step
171 
172  // here should be interpolation of the velocity at least if the interpolation time
173  // is not close to the solved_time of the water module
174  // for simplicity we use only last velocity field
175  if (pdata.velocity_changed) {
176  //DBGMSG("velocity update\n");
177  pdata.process->set_velocity_field( water->get_mh_dofhandler() );
178  pdata.velocity_changed = false;
179  }
180  if (pdata.process->time().tlevel() == 0) pdata.process->zero_time_step();
181 
182  pdata.process->update_solution();
183  pdata.process->output_data();
184  }
185 }
186 
187 
188 /**
189  * TODO:
190  * - have support for steady problems in TimeGovernor, make Noting problems steady
191  * - apply splitting of compute_one_step to particular models
192  * - how to set output time marks for steady problems (we need solved time == infinity) but
193  * add no time marks
194  * - allow create steady time governor without time marks (at least in nothing models)
195  * - pass refference to time marks in time governor constructor?
196  */
197 
199 {
200  START_TIMER("HC run simulation");
201  // following should be specified in constructor:
202  // value for velocity interpolation :
203  // theta = 0 velocity from beginning of transport interval (fully explicit method)
204  // theta = 0.5 velocity from center of transport interval ( mimic Crank-Nicholson)
205  // theta = 1.0 velocity from end of transport interval (partialy explicit scheme)
206  const double theta=0.5;
207 
208  {
209  START_TIMER("HC water zero time step");
210  water->zero_time_step();
211  }
212 
213 
214  // following cycle is designed to support independent time stepping of
215  // both processes. The question is which value of the water field use to compute a transport step.
216  // Meaningful cases are
217  // 1) beginning (fully explicit method)
218  // 2) center ( mimic Crank-Nicholson)
219  // 3) end of the interval (partialy explicit scheme)
220  // However with current implementation of the explicit transport on have to assembly transport matrix for
221  // every new value of the velocity field. So we have to keep same velocity field over some time interval t_dt
222  // which is further split into shorter time intervals ts_dt dictated by the CFL condition.
223  // One can consider t_dt as the transport time step and apply one of the previous three cases.
224  //
225  // The question is how to choose intervals t_dt. That should depend on variability of the velocity field in time.
226  // Currently we simply use t_dt == w_dt.
227 
228  is_end_all_=false;
229  while (! is_end_all_) {
230  is_end_all_ = true;
231 
232  double water_dt=water->time().estimate_dt();
233  if (water->time().is_end()) water_dt = TimeGovernor::inf_time;
234 
235  // in future here could be re-estimation of transport planed time according to
236  // evolution of the velocity field. Consider the case w_dt << t_dt and velocity almost constant in time
237  // which suddenly rise in time 3*w_dt. First we the planed transport time step t_dt could be quite big, but
238  // in time 3*w_dt we can reconsider value of t_dt to better capture changing velocity.
240  for(auto &pdata : processes_) {
241  pdata.process->set_time_upper_constraint(water_dt, "Flow time step");
242  pdata.velocity_time = theta * pdata.process->planned_time() + (1-theta) * pdata.process->solved_time();
243  min_velocity_time = min(min_velocity_time, pdata.velocity_time);
244  }
245 
246  // printing water and transport times every step
247  //MessageOut().fmt("HC_EXPL_SEQ: velocity_interpolation_time: {}, water_time: {} transport time: {}\n",
248  // velocity_interpolation_time, water->time().t(), transport_reaction->time().t());
249 
250  // if transport is off, transport should return infinity solved and planned times so that
251  // only water branch takes the place
252 
253  if (! water->time().is_end() ) {
254  is_end_all_=false;
255  if (water->solved_time() < min_velocity_time) {
256 
257  // solve water over the nearest transport interval
258  water->update_solution();
259 
260  // here possibly save solution from water for interpolation in time
261 
262  //water->time().view("WATER"); //show water time governor
263 
264  //water->output_data();
265  water->choose_next_time();
266 
267  for(auto &process : processes_) process.velocity_changed = true;
268 
269  }
270  }
271  advection_process_step(processes_[0]); // solute
272  advection_process_step(processes_[1]); // heat
273  }
274  //MessageOut().fmt("End of simulation at time: {}\n", max(solute->solved_time(), heat->solved_time()));
275 }
276 
277 
279  water.reset();
280  for(auto &pdata : processes_) pdata.process.reset();
281  delete mesh;
282 }
283 
284 
285 
286 
287 //-----------------------------------------------------------------------------
288 // vim: set cindent:
289 //-----------------------------------------------------------------------------
Mesh * mesh
mesh common to darcy flow and transport
Common abstract parent of all Field<...> classes.
Definition: field_common.hh:60
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
static const Input::Type::Record & get_input_type()
Definition: mesh.cc:54
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
Definition: mesh.h:97
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.
unsigned int n_elements() const
Definition: mesh.h:141
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:490
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 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:236
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:91
Class for declaration of the input data that are in string format.
Definition: type_base.hh:588
static const double inf_time
Infinity time used for steady case.
void set_task_info(string description, int size)