Flow123d  JS_before_hm-1013-g06f2edc
hm_iterative.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 hm_iterative.cc
15  * @brief
16  * @author Jan Stebel
17  */
18 
19 #include "hm_iterative.hh"
20 #include "system/sys_profiler.hh"
21 #include "input/input_type.hh"
22 #include "flow/richards_lmh.hh"
23 #include "fields/field_fe.hh" // for create_field_fe()
24 
25 
26 FLOW123D_FORCE_LINK_IN_CHILD(coupling_iterative)
27 
28 
29 namespace it = Input::Type;
30 
31 
32 const it::Record & HM_Iterative::get_input_type() {
33  return it::Record("Coupling_Iterative",
34  "Record with data for iterative coupling of flow and mechanics.\n")
38  .declare_key("flow_equation", RichardsLMH::get_input_type(),
40  "Flow equation, provides the velocity field as a result.")
41  .declare_key("mechanics_equation", Elasticity::get_input_type(),
42  "Mechanics, provides the displacement field.")
43  .declare_key("input_fields", it::Array(
45  .make_field_descriptor_type("Coupling_Iterative")),
47  "Input fields of the HM coupling.")
48  .declare_key( "iteration_parameter", it::Double(), it::Default("1"),
49  "Tuning parameter for iterative splitting. Its default value"
50  "corresponds to a theoretically optimal value with fastest convergence." )
51  .declare_key( "max_it", it::Integer(0), it::Default("100"),
52  "Maximal count of HM iterations." )
53  .declare_key( "min_it", it::Integer(0), it::Default("1"),
54  "Minimal count of HM iterations." )
55  .declare_key( "a_tol", it::Double(0), it::Default("0"),
56  "Absolute tolerance for difference in HM iteration." )
57  .declare_key( "r_tol", it::Double(0), it::Default("1e-7"),
58  "Relative tolerance for difference in HM iteration." )
59  .close();
60 }
61 
62 
63 const int HM_Iterative::registrar = Input::register_class< HM_Iterative, Mesh &, const Input::Record >("Coupling_Iterative")
65 
66 
68 {
69  *this += alpha.name("biot_alpha")
70  .description("Biot poroelastic coefficient.")
71  .units(UnitSI().dimensionless())
72  .input_default("0.0")
74 
75  *this += density.name("fluid_density")
76  .description("Volumetric mass density of the fluid.")
77  .units(UnitSI().kg().m(-3))
78  .input_default("0.0")
80 
81  *this += gravity.name("gravity")
82  .description("Gravitational acceleration constant.")
83  .units(UnitSI().m().s(-2))
84  .input_default("9.81")
86 
87  *this += beta.name("relaxation_beta")
88  .units(UnitSI().dimensionless())
90 
91  *this += pressure_potential.name("pressure_potential")
92  .units(UnitSI().m())
94 
95  *this += flow_source.name("extra_flow_source")
96  .units(UnitSI().s(-1))
98 }
99 
100 
102 {
103  // initialize coupling fields with FieldFE
104  set_mesh(mesh);
105 
106  potential_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(mesh, MixedPtr<FE_CR>());
107  pressure_potential.set_field(mesh.region_db().get_region_set("ALL"), potential_ptr_);
108 
109  beta_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(mesh, MixedPtr<FE_P_disc>(0));
110  beta.set_field(mesh.region_db().get_region_set("ALL"), beta_ptr_);
111 
112  flow_source_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(beta_ptr_->get_dofhandler());
113  flow_source.set_field(mesh.region_db().get_region_set("ALL"), flow_source_ptr_);
114 
115  old_pressure_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(beta_ptr_->get_dofhandler());
116  old_iter_pressure_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(beta_ptr_->get_dofhandler());
117  div_u_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(beta_ptr_->get_dofhandler());
118  old_div_u_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(beta_ptr_->get_dofhandler());
119 }
120 
121 
122 
124 : DarcyFlowInterface(mesh, in_record),
125  IterativeCoupling(in_record)
126 {
127  START_TIMER("HM constructor");
128  using namespace Input;
129 
130  time_ = new TimeGovernor(in_record.val<Record>("time"));
131  ASSERT( time_->is_default() == false ).error("Missing key 'time' in Coupling_Iterative.");
132 
133  // setup flow equation
134  Record flow_rec = in_record.val<Record>("flow_equation");
135  // Need explicit template types here, since reference is used (automatically passing by value)
136  flow_ = std::make_shared<RichardsLMH>(*mesh_, flow_rec, time_);
137  flow_->initialize();
138  std::stringstream ss; // print warning message with table of uninitialized fields
139  if ( FieldCommon::print_message_table(ss, "flow") )
140  WarningOut() << ss.str();
141 
142  // setup mechanics
143  Record mech_rec = in_record.val<Record>("mechanics_equation");
144  mechanics_ = std::make_shared<Elasticity>(*mesh_, mech_rec, this->time_);
145  mechanics_->data()["cross_section"].copy_from(flow_->data()["cross_section"]);
146  mechanics_->initialize();
147 
148  // read parameters controlling the iteration
149  beta_ = in_record.val<double>("iteration_parameter");
150 
151  this->eq_data_ = &data_;
152 
153  // setup input fields
154  data_.set_input_list( in_record.val<Input::Array>("input_fields"), time() );
155 
156  data_.initialize(*mesh_);
157  mechanics_->set_potential_load(data_.pressure_potential);
158 }
159 
160 
162 {
163 }
164 
165 
166 template<int dim, class Value>
167 void copy_field(const FieldCommon &from_field_common, FieldFE<dim, Value> &to_field)
168 {
169  auto dh = to_field.get_dofhandler();
170  auto vec = to_field.vec();
171  Field<dim,Value> from_field;
172  from_field.copy_from(from_field_common);
173 
174  for ( auto cell : dh->own_range() )
175  vec[cell.local_idx()] = from_field.value(cell.elm().centre(), cell.elm());
176 }
177 
178 
179 
181 {
183  std::stringstream ss;
184  if ( FieldCommon::print_message_table(ss, "coupling_iterative") )
185  WarningOut() << ss.str();
186 
187  flow_->zero_time_step();
189  mechanics_->zero_time_step();
190 
191  copy_field(*flow_->data().field("pressure_p0"), *data_.old_pressure_ptr_);
192  copy_field(*flow_->data().field("pressure_p0"), *data_.old_iter_pressure_ptr_);
193  copy_field(mechanics_->data().output_divergence, *data_.old_div_u_ptr_);
194  copy_field(mechanics_->data().output_divergence, *data_.div_u_ptr_);
195 }
196 
197 
199 {
200  time_->next_time();
201  time_->view("HM");
203 
204  solve_step();
205 }
206 
208 {
209  // pass displacement (divergence) to flow
210  // and solve flow problem
212  flow_->solve_time_step(false);
213 
214  // pass pressure to mechanics and solve mechanics
216  mechanics_->solve_linear_system();
217 }
218 
219 
221 {
222  mechanics_->update_output_fields();
223  copy_field(mechanics_->data().output_divergence, *data_.div_u_ptr_);
224  copy_field(*flow_->data().field("pressure_p0"), *data_.old_iter_pressure_ptr_);
225 }
226 
227 
229 {
230  flow_->accept_time_step();
231  flow_->output_data();
232  mechanics_->output_data();
233 
234  copy_field(*flow_->data().field("pressure_p0"), *data_.old_pressure_ptr_);
235  copy_field(mechanics_->data().output_divergence, *data_.old_div_u_ptr_);
236 }
237 
238 
240 {
241  auto potential_vec_ = data_.potential_ptr_->vec();
242  auto dh = data_.potential_ptr_->get_dofhandler();
243  Field<3, FieldValue<3>::Scalar> field_edge_pressure;
244  field_edge_pressure.copy_from(*flow_->data().field("pressure_edge"));
245  for ( auto ele : dh->local_range() )
246  {
247  auto elm = ele.elm();
248  LocDofVec dof_indices = ele.get_loc_dof_indices();
249  for ( auto side : ele.side_range() )
250  {
251  double alpha = data_.alpha.value(side.centre(), elm);
252  double density = data_.density.value(side.centre(), elm);
253  double gravity = data_.gravity.value(side.centre(), elm);
254  double pressure = field_edge_pressure.value(side.centre(), elm);
255  double potential = -alpha*density*gravity*pressure;
256 
257  potential_vec_[dof_indices[side.side_idx()]] = potential;
258  }
259  }
260 
262  mechanics_->set_potential_load(data_.pressure_potential);
263 }
264 
265 
267 {
268  auto beta_vec = data_.beta_ptr_->vec();
269  auto src_vec = data_.flow_source_ptr_->vec();
270  auto dh = data_.beta_ptr_->get_dofhandler();
271  Field<3,FieldValue<3>::Scalar> field_ele_pressure;
272  field_ele_pressure.copy_from(*flow_->data().field("pressure_p0"));
273  for ( auto ele : dh->local_range() )
274  {
275  auto elm = ele.elm();
276 
277  double alpha = data_.alpha.value(elm.centre(), elm);
278  double young = mechanics_->data().young_modulus.value(elm.centre(), elm);
279  double poisson = mechanics_->data().poisson_ratio.value(elm.centre(), elm);
280  double beta = beta_ * 0.5*alpha*alpha/(2*lame_mu(young, poisson)/elm.dim() + lame_lambda(young, poisson));
281 
282  double old_p = data_.old_pressure_ptr_->value(elm.centre(), elm);
283  double p = field_ele_pressure.value(elm.centre(), elm);
284  double div_u = data_.div_u_ptr_->value(elm.centre(), elm);
285  double old_div_u = data_.old_div_u_ptr_->value(elm.centre(), elm);
286  double src = (beta*(p-old_p) + alpha*(old_div_u - div_u)) / time_->dt();
287 
288  beta_vec[ele.local_idx()] = beta;
289  src_vec[ele.local_idx()] = src;
290  }
291 
294  flow_->set_extra_storativity(data_.beta);
295  flow_->set_extra_source(data_.flow_source);
296 }
297 
298 
299 void HM_Iterative::compute_iteration_error(double& abs_error, double& rel_error)
300 {
301  auto dh = data_.beta_ptr_->get_dofhandler();
302  double p_dif2 = 0, p_norm2 = 0;
303  Field<3,FieldValue<3>::Scalar> field_ele_pressure;
304  field_ele_pressure.copy_from(*flow_->data().field("pressure_p0"));
305  for (auto cell : dh->own_range())
306  {
307  auto elm = cell.elm();
308  double new_p = field_ele_pressure.value(elm.centre(), elm);
309  double old_p = data_.old_iter_pressure_ptr_->value(elm.centre(), elm);
310  p_dif2 += pow(new_p - old_p, 2)*elm.measure();
311  p_norm2 += pow(old_p, 2)*elm.measure();
312  }
313 
314  double send_data[] = { p_dif2, p_norm2 };
315  double recv_data[2];
316  MPI_Allreduce(&send_data, &recv_data, 2, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
317  abs_error = sqrt(recv_data[0]);
318  rel_error = abs_error / sqrt(recv_data[1]);
319 
320  MessageOut().fmt("HM Iteration {} abs. difference: {} rel. difference: {}\n"
321  "--------------------------------------------------------",
322  iteration(), abs_error, rel_error);
323 
324  if(iteration() >= max_it_ && (abs_error > a_tol_ || rel_error > r_tol_))
325  MessageOut().fmt("HM solver did not converge in {} iterations.\n", iteration());
326 }
327 
328 
329 
331  flow_.reset();
332  mechanics_.reset();
333 }
334 
335 
336 
TimeGovernor & time()
Definition: equation.hh:149
FieldSet * eq_data_
Definition: equation.hh:227
Common abstract parent of all Field<...> classes.
Definition: field_common.hh:73
Accessor to input data conforming to declared Array.
Definition: accessors.hh:566
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
ArmaVec< double, N > vec
Definition: armor.hh:861
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > old_div_u_ptr_
Class for solution of fully coupled flow and mechanics using fixed-stress iterative splitting...
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:602
double beta_
Tuning parameter for iterative splitting.
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
static const int registrar
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
#define MessageOut()
Macro defining &#39;message&#39; record of log.
Definition: logger.hh:255
void update_potential()
Abstract linear system class.
Definition: balance.hh:40
void compute_iteration_error(double &abs_error, double &rel_error) override
Compute absolute and relative error in the solution.
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:92
void initialize(Mesh &mesh)
unsigned int iteration()
Definition: hm_iterative.hh:78
void next_time()
Proceed to the next time according to current estimated time step.
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
Definition: elasticity.cc:43
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
double a_tol_
Absolute tolerance for difference between two succeeding iterations.
#define MPI_SUM
Definition: mpi.h:196
Definition: mesh.h:78
void update_flow_fields()
double lame_mu(double young, double poisson)
Definition: elasticity.cc:128
void zero_time_step() override
static const Input::Type::Record & record_template()
Definition: hm_iterative.hh:39
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > potential_ptr_
FieldFE for pressure_potential field.
const TimeStep & step(int index=-1) const
Class for declaration of the integral input data.
Definition: type_base.hh:483
Basic time management functionality for unsteady (and steady) solvers (class Equation).
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
static const Input::Type::Record & get_input_type()
Define input record.
Definition: hm_iterative.cc:32
void update_after_iteration() override
Save data (e.g. solution fields) for the next iteration.
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
Class for declaration of inputs sequences.
Definition: type_base.hh:339
void view(const char *name="") const
static Input::Type::Record & record_template()
Template Record with common keys for derived equations.
Definition: equation.cc:35
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > old_pressure_ptr_
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > div_u_ptr_
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
void update_solution() override
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:534
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
Mesh & mesh()
Definition: equation.hh:177
Field< 3, FieldValue< 3 >::Scalar > pressure_potential
Potential -alpha*pressure whose gradient is passed to mechanics as additional load.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > flow_source_ptr_
FieldCommon & input_default(const string &input_default)
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Definition: field_flag.hh:58
void set_time_result_changed()
Manually mark flag that the field has been changed.
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
Field< 3, FieldValue< 3 >::Scalar > alpha
Biot coefficient.
const Ret val(const string &key) const
std::shared_ptr< Elasticity > mechanics_
solute transport with chemistry through operator splitting
#define START_TIMER(tag)
Starts a timer with specified tag.
Mesh * mesh_
Definition: equation.hh:218
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:434
Field< 3, FieldValue< 3 >::Scalar > gravity
Standard gravity.
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
static Input::Type::Abstract & get_input_type()
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh:51
Field< 3, FieldValue< 3 >::Scalar > flow_source
void initialize() override
Field< 3, FieldValue< 3 >::Scalar > beta
#define MPI_DOUBLE
Definition: mpi.h:156
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
std::shared_ptr< DOFHandlerMultiDim > get_dofhandler() const
Definition: field_fe.hh:142
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
Definition: mpi.h:612
FieldCommon & description(const string &description)
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
Definition: field_set.hh:193
static const Input::Type::Record & get_input_type()
Definition: richards_lmh.cc:79
double dt() const
void set_field(const RegionSet &domain, FieldBasePtr field, double time=0.0)
Definition: field.impl.hh:227
unsigned int max_it_
Maximal number of iterations.
HM_Iterative(Mesh &mesh, Input::Record in_record)
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:157
void solve_iteration() override
Solve equations and update data (fields).
double lame_lambda(double young, double poisson)
Definition: elasticity.cc:134
void copy_field(const FieldCommon &from_field_common, FieldFE< dim, Value > &to_field)
#define WarningOut()
Macro defining &#39;warning&#39; record of log.
Definition: logger.hh:258
FieldCommon & name(const string &name)
VectorMPI & vec()
Definition: field_fe.hh:146
void update_after_converged() override
Save data after iterations have finished.
double r_tol_
Relative tolerance for difference between two succeeding iterations.
void set_mesh(const Mesh &mesh)
Definition: field_set.hh:186
Record type proxy class.
Definition: type_record.hh:182
FieldCommon & flags(FieldFlag::Flags::Mask mask)
void copy_from(const FieldCommon &other) override
Definition: field.impl.hh:345
Class for representation SI units of Fields.
Definition: unit_si.hh:40
std::shared_ptr< RichardsLMH > flow_
steady or unsteady water flow simulator based on MH scheme
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:96
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > old_iter_pressure_ptr_
Field< 3, FieldValue< 3 >::Scalar > density
Density of fluid.
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
Definition: global_defs.h:180
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > beta_ptr_
TimeGovernor * time_
Definition: equation.hh:219
Definition: field.hh:60