Flow123d  jenkins-Flow123d-windows32-release-multijob-51
hc_explicit_sequential.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief
27  *
28  * @author Jan Brezina
29  */
30 
32 #include "flow/darcy_flow_mh.hh"
35 #include "transport/transport.h"
37 #include "mesh/mesh.h"
38 #include "mesh/msh_gmshreader.h"
39 #include "system/sys_profiler.hh"
40 #include "input/input_type.hh"
42 #include "transport/heat_model.hh"
43 
44 
45 namespace it = Input::Type;
46 
48  = it::AbstractRecord("Problem",
49  "The root record of description of particular the problem to solve.")
50  .declare_key("description",it::String(),
51  "Short description of the solved problem.\n"
52  "Is displayed in the main log, and possibly in other text output files.")
54  "Computational mesh common to all equations.");
55 
56 
58  = it::Record("SequentialCoupling",
59  "Record with data for a general sequential coupling.\n")
60  .derive_from( CouplingBase::input_type )
62  "Simulation time frame and time step.")
64  "Primary equation, have all data given.")
65  .declare_key("secondary_equation", TransportBase::input_type,
66  "The equation that depends (the velocity field) on the result of the primary equation.");
67 
68 
69 
70 
71 /**
72  * FUNCTION "MAIN" FOR COMPUTING MIXED-HYBRID PROBLEM FOR UNSTEADY SATURATED FLOW
73  */
75 {
76  START_TIMER("HC constructor");
77  using namespace Input;
78 
79  // Read mesh
80  {
81  mesh = new Mesh( in_record.val<Record>("mesh") );
83 
84  //getting description for the Profiler
85  string description;
86  in_record.opt_val<string>("description", description);
87 
89  description,
90  //"Description has to be set in main. by different method.",
91  mesh->n_elements());
92  }
93 
94  // setup primary equation - water flow object
95  AbstractRecord prim_eq = in_record.val<AbstractRecord>("primary_equation");
96  if (prim_eq.type() == DarcyFlowMH_Steady::input_type ) {
97  water = new DarcyFlowMH_Steady(*mesh, prim_eq);
98  } else if (prim_eq.type() == DarcyFlowMH_Unsteady::input_type ) {
99  water = new DarcyFlowMH_Unsteady(*mesh, prim_eq);
100  } else if (prim_eq.type() == DarcyFlowLMH_Unsteady::input_type ) {
101  water = new DarcyFlowLMH_Unsteady(*mesh, prim_eq);
102  } else {
103  xprintf(UsrErr,"Equation type not implemented.");
104  }
105 
106 
107 
108  // TODO: optionally setup transport objects
109  Iterator<AbstractRecord> it = in_record.find<AbstractRecord>("secondary_equation");
110  if (it) {
111  if (it->type() == TransportOperatorSplitting::input_type)
112  {
114  }
116  {
118  }
119  else if (it->type() == TransportDG<HeatTransferModel>::input_type)
120  {
122  }
123  else
124  {
125  xprintf(PrgErr,"Value of TYPE in the Transport an AbstractRecord out of set of descendants.\n");
126  }
127 
128  // setup fields
129  transport_reaction->data()["cross_section"]
130  .copy_from(water->data()["cross_section"]);
131 
132  } else {
134  }
135 }
136 
137 
138 
139 /**
140  * TODO:
141  * - have support for steady problems in TimeGovernor, make Noting problems steady
142  * - apply splitting of compute_one_step to particular models
143  * - how to set output time marks for steady problems (we need solved time == infinity) but
144  * add no time marks
145  * - allow create steady time governor without time marks (at least in nothing models)
146  * - pass refference to time marks in time governor constructor?
147  */
148 
150 {
151  START_TIMER("HC run simulation");
152  // following should be specified in constructor:
153  // value for velocity interpolation :
154  // theta = 0 velocity from beginning of transport interval (fully explicit method)
155  // theta = 0.5 velocity from center of transport interval ( mimic Crank-Nicholson)
156  // theta = 1.0 velocity from end of transport interval (partialy explicit scheme)
157  const double theta=0.5;
158 
159 
160  double velocity_interpolation_time;
161  bool velocity_changed=true;
162 
163 
165 
166 
167 
168  // following cycle is designed to support independent time stepping of
169  // both processes. The question is which value of the water field use to compute a transport step.
170  // Meaningful cases are
171  // 1) beginning (fully explicit method)
172  // 2) center ( mimic Crank-Nicholson)
173  // 3) end of the interval (partialy explicit scheme)
174  // However with current implementation of the explicit transport on have to assembly transport matrix for
175  // every new value of the velocity field. So we have to keep same velocity field over some time interval t_dt
176  // which is further split into shorter time intervals ts_dt dictated by the CFL condition.
177  // One can consider t_dt as the transport time step and apply one of the previous three cases.
178  //
179  // The question is how to choose intervals t_dt. That should depend on variability of the velocity field in time.
180  // Currently we simply use t_dt == w_dt.
181 
182  while (! (water->time().is_end() && transport_reaction->time().is_end() ) ) {
183 
185  // in future here could be re-estimation of transport planed time according to
186  // evolution of the velocity field. Consider the case w_dt << t_dt and velocity almost constant in time
187  // which suddenly rise in time 3*w_dt. First we the planed transport time step t_dt could be quite big, but
188  // in time 3*w_dt we can reconsider value of t_dt to better capture changing velocity.
189  velocity_interpolation_time= theta * transport_reaction->planned_time() + (1-theta) * transport_reaction->solved_time();
190 
191  // printing water and transport times every step
192  //xprintf(Msg,"HC_EXPL_SEQ: velocity_interpolation_time: %f, water_time: %f transport time: %f\n",
193  // velocity_interpolation_time, water->time().t(), transport_reaction->time().t());
194 
195  // if transport is off, transport should return infinity solved and planned times so that
196  // only water branch takes the place
197  if (water->solved_time() < velocity_interpolation_time) {
198  // solve water over the nearest transport interval
200 
201  // here possibly save solution from water for interpolation in time
202 
203  //water->time().view("WATER"); //show water time governor
204 
205  //water->output_data();
207 
208  velocity_changed = true;
209  } else {
210  // having information about velocity field we can perform transport step
211 
212  // here should be interpolation of the velocity at least if the interpolation time
213  // is not close to the solved_time of the water module
214  // for simplicity we use only last velocity field
215  if (velocity_changed) {
216  //DBGMSG("velocity update\n");
218  velocity_changed = false;
219  }
221 
223 
224  //transport_reaction->time().view("TRANSPORT"); //show transport time governor
225 
227  }
228 
229  }
230  xprintf(Msg, "End of simulation at time: %f\n", transport_reaction->solved_time());
231 }
232 
233 
235  delete water;
236  delete transport_reaction;
237  delete mesh;
238 }
239 
240 
241 
242 //-----------------------------------------------------------------------------
243 // vim: set cindent:
244 //-----------------------------------------------------------------------------
FieldSet & data()
Definition: equation.hh:188
virtual void zero_time_step()
Definition: equation.hh:97
Mesh * mesh
mesh common to darcy flow and transport
Output class for darcy_flow_mh model.
static Input::Type::AbstractRecord input_type
Common specification of the input record for secondary equations.
Definition: system.hh:72
int tlevel() const
Transport with dispersion implemented using discontinuous Galerkin method.
static Profiler * instance()
bool is_end() const
Returns true if the actual time is greater than or equal to the end time.
static Input::Type::AbstractRecord input_type
DarcyFlowMH * water
steady or unsteady water flow simulator based on MH scheme
static Default obligatory()
Definition: type_record.hh:87
???
Coupling of a transport model with a reaction model by operator splitting.
Definition: mesh.h:108
Iterator< Ret > find(const string &key) const
static Input::Type::Record input_type
Record & derive_from(AbstractRecord &parent)
Definition: type_record.cc:169
const MH_DofHandler & get_mh_dofhandler()
static Input::Type::Record input_type
Discontinuous Galerkin method for equation of transport with dispersion.
virtual void choose_next_time()
Definition: equation.hh:125
virtual void output_data()
Write computed fields.
Definition: equation.hh:211
AdvectionProcessBase * transport_reaction
explicit transport with chemistry through operator splitting
HC_ExplicitSequential(Input::Record in_record)
static Default optional()
Definition: type_record.hh:100
Mixed-hybrid solution of unsteady Darcy flow.
bool opt_val(const string &key, Ret &value) const
unsigned int n_elements() const
Definition: mesh.h:137
Accessor to the data with type Type::Record.
Definition: accessors.hh:308
const Ret val(const string &key) const
#define xprintf(...)
Definition: system.hh:100
#define START_TIMER(tag)
Starts a timer with specified tag.
Class for declaration of polymorphic Record.
Definition: type_record.hh:463
static Input::Type::Record input_type
AbstractRecord & declare_key(const string &key, const KeyType &type, const Default &default_value, const string &description)
Definition: type_record.cc:422
static Input::Type::Record input_type
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:423
Discontinuous Galerkin method for equation of transport with dispersion.
Empty transport class.
Definition: system.hh:72
static Input::Type::AbstractRecord input_type
virtual void update_solution()
Definition: equation.hh:105
static Input::Type::Record input_type
Definition: mesh.h:111
double dt() const
Definition: system.hh:72
TimeGovernor const & time()
Definition: equation.hh:143
Edge lumped mixed-hybrid solution of unsteady Darcy flow.
double solved_time()
Definition: equation.hh:165
Mixed-hybrid of steady Darcy flow with sources and variable density.
virtual void set_velocity_field(const MH_DofHandler &dh)=0
Discontinuous Galerkin method for equation of transport with dispersion.
static Input::Type::Record input_type
mixed-hybrid model of linear Darcy flow, possibly unsteady.
Record type proxy class.
Definition: type_record.hh:161
static Input::Type::Record input_type
Declare input record type for the equation TransportOperatorSplittiong.
virtual void set_time_upper_constraint(double dt)
Definition: equation.hh:131
void init_from_input()
Definition: mesh.cc:198
void set_task_info(string description, int size)
Record & declare_key(const string &key, const KeyType &type, const Default &default_value, const string &description)
Definition: type_record.cc:386
double planned_time()
Definition: equation.hh:159