Flow123d  3.9.0-3aaaea010
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/darcy_flow_lmh.hh"
23 #include "fields/field_fe.hh" // for create_field_fe()
24 #include "fields/field_model.hh" // for Model
25 #include "assembly_hm.hh"
26 
27 
28 FLOW123D_FORCE_LINK_IN_CHILD(coupling_iterative)
29 
30 
31 namespace it = Input::Type;
32 
33 
35  inline double operator() (double alpha, double density, double gravity, double pressure)
36  {
37  return -alpha*density*gravity*pressure;
38  }
39 };
40 
41 
43 
44  fn_hm_coupling_beta(double beta_f) : beta_factor(beta_f) {}
45 
46 
47  inline double operator() (double alpha, double lame_mu, double lame_lambda, double density, double gravity)
48  {
49  return beta_factor*0.5*alpha*alpha/(2*lame_mu/3 + lame_lambda)*density*gravity;
50  }
51 
52 private:
53 
54  const double beta_factor; ///< User-defined factor for iteration parameter.
55 
56 };
57 
58 
60 
61  fn_fluid_source(const TimeGovernor *time) : time_(time) {}
62 
63  inline double operator() (double alpha, double beta, double pressure, double old_pressure, double div_u, double old_div_u)
64  {
65  return (beta*(pressure-old_pressure) + alpha*(old_div_u - div_u)) / time_->dt();
66  }
67 
68 private:
69 
71 };
72 
73 
75  return it::Record("Coupling_Iterative",
76  "Record with data for iterative coupling of flow and mechanics.\n")
80  .declare_key("flow_equation", DarcyLMH::get_input_type(),
82  "Flow equation, provides the velocity field as a result.")
83  .declare_key("mechanics_equation", Elasticity::get_input_type(),
85  "Mechanics, provides the displacement field.")
86  .declare_key("input_fields", it::Array(
88  .make_field_descriptor_type("Coupling_Iterative")),
90  "Input fields of the HM coupling.")
91  .declare_key( "iteration_parameter", it::Double(), it::Default("1"),
92  "Tuning parameter for iterative splitting. Its default value "
93  "corresponds to a theoretically optimal value with fastest convergence." )
94  .declare_key( "max_it", it::Integer(0), it::Default("100"),
95  "Maximal count of HM iterations." )
96  .declare_key( "min_it", it::Integer(0), it::Default("1"),
97  "Minimal count of HM iterations." )
98  .declare_key( "a_tol", it::Double(0), it::Default("0"),
99  "Absolute tolerance for difference in HM iteration." )
100  .declare_key( "r_tol", it::Double(0), it::Default("1e-7"),
101  "Relative tolerance for difference in HM iteration." )
102  .close();
103 }
104 
105 
106 const int HM_Iterative::registrar = Input::register_class< HM_Iterative, Mesh &, const Input::Record >("Coupling_Iterative")
108 
109 
111 {
112  *this += alpha.name("biot_alpha")
113  .description("Biot poroelastic coefficient.")
114  .units(UnitSI().dimensionless())
115  .input_default("0.0")
117 
118  *this += density.name("fluid_density")
119  .description("Volumetric mass density of the fluid.")
120  .units(UnitSI().kg().m(-3))
121  .input_default("0.0")
123 
124  *this += gravity.name("gravity")
125  .description("Gravitational acceleration constant.")
126  .units(UnitSI().m().s(-2))
127  .input_default("9.81")
129 
130  *this += beta.name("relaxation_beta")
131  .description("Parameter of numerical method for iterative solution of hydro-mechanical coupling.")
132  .units(UnitSI().dimensionless())
134 
135  *this += pressure_potential.name("pressure_potential")
136  .description("Coupling term entering the mechanics equation.")
137  .units(UnitSI().m())
139 
140  *this += old_pressure.name("old_pressure")
141  .description("Pressure from last computed time.")
142  .units(UnitSI().m())
144 
145  *this += old_iter_pressure.name("old_iter_pressure")
146  .description("Pressure from last computed iteration.")
147  .units(UnitSI().m())
149 
150  *this += old_div_u.name("old_displacement_divergence")
151  .description("Displacement divergence from last computed time.")
152  .units(UnitSI().dimensionless())
154 
155  *this += ref_pressure_potential.name("ref_pressure_potential")
156  .description("Pressure potential on boundary (taking into account the flow boundary condition.")
157  .units(UnitSI().m())
159 
160  *this += flow_source.name("extra_flow_source")
161  .description("Coupling term entering the flow equation.")
162  .units(UnitSI().s(-1))
164 }
165 
166 
168 {
169  // initialize coupling fields with FieldFE
170  set_mesh(mesh);
171 
172  pressure_potential.set(Model<3, FieldValue<3>::Scalar>::create(
174  alpha,
175  density,
176  gravity,
177  eq_data.flow_->eq_fields().field_edge_pressure
178  ), 0.0);
179 
180  ref_potential_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(mesh, MixedPtr<FE_CR>());
181  ref_pressure_potential.set(ref_potential_ptr_, 0.0);
182 
183  beta.set(Model<3, FieldValue<3>::Scalar>::create(
184  fn_hm_coupling_beta(beta_),
185  alpha,
186  eq_data.mechanics_->eq_fields().lame_mu,
187  eq_data.mechanics_->eq_fields().lame_lambda,
188  density,
189  gravity
190  ), 0.0);
191 
192  auto old_pressure_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(eq_data.flow_->eq_data().dh_cr_, &eq_data.flow_->eq_data().p_edge_solution_previous_time);
193  old_pressure.set(old_pressure_ptr_, 0.0);
194 
195  old_iter_pressure_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(mesh, MixedPtr<FE_P_disc>(0));
196  old_iter_pressure.set(old_iter_pressure_ptr_, 0.0);
197 
198  old_div_u_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(eq_data.mechanics_->eq_fields().output_div_ptr->get_dofhandler());
199  old_div_u.set(old_div_u_ptr_, 0.0);
200 
201  flow_source.set(Model<3, FieldValue<3>::Scalar>::create(
203  alpha,
204  beta,
205  eq_data.flow_->eq_fields().field_edge_pressure,
206  old_pressure,
207  eq_data.mechanics_->eq_fields().output_divergence,
208  old_div_u
209  ), 0.0);
210 }
211 
212 
213 
215 : DarcyFlowInterface(mesh, in_record),
216  IterativeCoupling(in_record),
217  flow_potential_assembly_(nullptr),
218  residual_assembly_(nullptr)
219 {
220  START_TIMER("HM constructor");
221  using namespace Input;
222 
223  time_ = new TimeGovernor(in_record.val<Record>("time"));
224  ASSERT( time_->is_default() == false ).error("Missing key 'time' in Coupling_Iterative.");
225 
226  // setup flow equation
227  Record flow_rec = in_record.val<Record>("flow_equation");
228  // Need explicit template types here, since reference is used (automatically passing by value)
229  eq_data_.flow_ = std::make_shared<DarcyLMH>(*mesh_, flow_rec, time_);
230 
231  // setup mechanics
232  Record mech_rec = in_record.val<Record>("mechanics_equation");
233  eq_data_.mechanics_ = std::make_shared<Elasticity>(*mesh_, mech_rec, this->time_);
234  eq_data_.mechanics_->initialize();
235 
236  // setup coupling fields and finish initialization of flow
237  eq_data_.mechanics_->eq_fields()["cross_section"].copy_from(eq_data_.flow_->eq_fields()["cross_section"]);
238  eq_data_.flow_->eq_fields() += eq_data_.mechanics_->eq_fields()["cross_section_updated"];
239  eq_data_.flow_->eq_fields() += eq_data_.mechanics_->eq_fields()["stress"];
240  eq_data_.flow_->eq_fields() += eq_data_.mechanics_->eq_fields()["von_mises_stress"];
241  eq_data_.flow_->eq_fields() += eq_data_.mechanics_->eq_fields()["mean_stress"];
242  eq_data_.flow_->initialize();
243  std::stringstream ss; // print warning message with table of uninitialized fields
244  if ( FieldCommon::print_message_table(ss, "flow") )
245  WarningOut() << ss.str();
246 
247  eq_fields_ += *eq_data_.flow_->eq_fields().field("pressure_edge");
248 
249  this->eq_fieldset_ = &eq_fields_;
250 
251  // setup input fields
252  eq_fields_.set_input_list( in_record.val<Input::Array>("input_fields"), time() );
253 
254  eq_fields_.initialize(*mesh_, eq_data_, time_, input_record_.val<double>("iteration_parameter"));
256 
258 }
259 
260 
262 {
265 }
266 
267 
268 template<int dim, class Value>
269 void copy_field(const FieldCommon &from_field_common, FieldFE<dim, Value> &to_field)
270 {
271  auto dh = to_field.get_dofhandler();
272  auto vec = to_field.vec();
273  Field<dim,Value> from_field;
274  from_field.copy_from(from_field_common);
275 
276  for ( auto cell : dh->own_range() )
277  vec.set( cell.local_idx(), from_field.value(cell.elm().centre(), cell.elm()) );
278 }
279 
280 
281 
283 {
285  std::stringstream ss;
286  if ( FieldCommon::print_message_table(ss, "coupling_iterative") )
287  WarningOut() << ss.str();
288 
289  eq_data_.mechanics_->update_output_fields(); // init field values for use in flow
290  eq_data_.flow_->zero_time_step();
292  eq_data_.mechanics_->zero_time_step();
293 
294  copy_field(*eq_data_.flow_->eq_fields().field("pressure_p0"), *eq_fields_.old_iter_pressure_ptr_);
295  copy_field(eq_data_.mechanics_->eq_fields().output_divergence, *eq_fields_.old_div_u_ptr_);
297 }
298 
299 
301 {
302  time_->next_time();
303  time_->view("HM");
305 
306  solve_step();
307 }
308 
310 {
311  // pass displacement (divergence) to flow
312  // and solve flow problem
314  eq_data_.flow_->solve_time_step(false);
315 
316  // pass pressure to mechanics and solve mechanics
318  eq_data_.mechanics_->solve_linear_system();
319 }
320 
321 
323 {
324  eq_data_.mechanics_->update_output_fields();
325  copy_field(*eq_data_.flow_->eq_fields().field("pressure_p0"), *eq_fields_.old_iter_pressure_ptr_);
327 }
328 
329 
331 {
332  eq_data_.flow_->accept_time_step();
333  eq_data_.flow_->output_data();
334  eq_data_.mechanics_->output_data();
335 
336  copy_field(eq_data_.mechanics_->eq_fields().output_divergence, *eq_fields_.old_div_u_ptr_);
337 }
338 
339 
341 {
342  auto ref_potential_vec_ = eq_fields_.ref_potential_ptr_->vec();
343  auto dh = eq_fields_.ref_potential_ptr_->get_dofhandler();
344 
345  ref_potential_vec_.zero_entries();
347 
348  ref_potential_vec_.local_to_ghost_begin();
349  ref_potential_vec_.local_to_ghost_end();
353 }
354 
355 
357 {
360  eq_data_.flow_->set_extra_storativity(eq_fields_.beta);
361  eq_data_.flow_->set_extra_source(eq_fields_.flow_source);
362 }
363 
364 
365 void HM_Iterative::compute_iteration_error(double& abs_error, double& rel_error)
366 {
367  auto dh = eq_fields_.old_iter_pressure_ptr_->get_dofhandler();
368  eq_data_.p_dif2 = 0;
369  eq_data_.p_norm2 = 0;
370 
372 
373  double send_data[] = { eq_data_.p_dif2, eq_data_.p_norm2 };
374  double recv_data[2];
375  MPI_Allreduce(&send_data, &recv_data, 2, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
376  abs_error = sqrt(recv_data[0]);
377  rel_error = abs_error / (sqrt(recv_data[1]) + std::numeric_limits<double>::min());
378 
379  MessageOut().fmt("HM Iteration {} abs. difference: {} rel. difference: {}\n"
380  "--------------------------------------------------------",
381  iteration(), abs_error, rel_error);
382 
383  if(iteration() >= max_it_ && (abs_error > a_tol_ || rel_error > r_tol_))
384  THROW(ExcSolverDiverge() << EI_Reason("Reached max_it."));
385 }
386 
387 
388 
390  eq_data_.flow_.reset();
391  eq_data_.mechanics_.reset();
392  if (flow_potential_assembly_ != nullptr) delete flow_potential_assembly_;
393  if (residual_assembly_ != nullptr) delete residual_assembly_;
394 }
395 
396 
397 
FieldCommon::units
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
Definition: field_common.hh:152
LimitSide::right
@ right
HM_Iterative::EqFields::alpha
Field< 3, FieldValue< 3 >::Scalar > alpha
Biot coefficient.
Definition: hm_iterative.hh:157
HM_Iterative::get_input_type
static const Input::Type::Record & get_input_type()
Define input record.
Definition: hm_iterative.cc:74
FieldFlag::in_rhs
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Definition: field_flag.hh:51
IterativeCoupling::solve_step
void solve_step()
Definition: hm_iterative.hh:72
EquationBase::mesh_
Mesh * mesh_
Definition: equation.hh:219
HM_Iterative::EqFields::EqFields
EqFields()
Definition: hm_iterative.cc:110
HM_Iterative::EqFields::old_div_u_ptr_
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > old_div_u_ptr_
Definition: hm_iterative.hh:173
IterativeCoupling::iteration
unsigned int iteration()
Definition: hm_iterative.hh:88
HM_Iterative::EqFields::beta
Field< 3, FieldValue< 3 >::Scalar > beta
Definition: hm_iterative.hh:160
HM_Iterative::eq_fields_
EqFields eq_fields_
Definition: hm_iterative.hh:204
fn_pressure_potential::operator()
double operator()(double alpha, double density, double gravity, double pressure)
Definition: hm_iterative.cc:35
Armor::vec
ArmaVec< double, N > vec
Definition: armor.hh:885
TimeGovernor::dt
double dt() const
Definition: time_governor.hh:565
HM_Iterative::EqFields::ref_potential_ptr_
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > ref_potential_ptr_
FieldFE for pressure_potential field.
Definition: hm_iterative.hh:171
HM_Iterative::EqFields::pressure_potential
Field< 3, FieldValue< 3 >::Scalar > pressure_potential
Potential -alpha*pressure whose gradient is passed to mechanics as additional load.
Definition: hm_iterative.hh:163
FieldSet::set_time
bool set_time(const TimeStep &time, LimitSide limit_side)
Definition: field_set.cc:244
fn_pressure_potential
Definition: hm_iterative.cc:34
fn_hm_coupling_beta::beta_factor
const double beta_factor
User-defined factor for iteration parameter.
Definition: hm_iterative.cc:54
HM_Iterative::EqFields::density
Field< 3, FieldValue< 3 >::Scalar > density
Density of fluid.
Definition: hm_iterative.hh:158
Input
Abstract linear system class.
Definition: balance.hh:40
HM_Iterative::EqFields::flow_source
Field< 3, FieldValue< 3 >::Scalar > flow_source
Definition: hm_iterative.hh:165
ASSERT
#define ASSERT(expr)
Definition: asserts.hh:351
EquationBase::time_
TimeGovernor * time_
Definition: equation.hh:220
Input::Type::Integer
Class for declaration of the integral input data.
Definition: type_base.hh:483
IterativeCoupling
Definition: hm_iterative.hh:42
Input::Record::val
const Ret val(const string &key) const
Definition: accessors_impl.hh:31
HM_Iterative::EqFields::old_iter_pressure
Field< 3, FieldValue< 3 >::Scalar > old_iter_pressure
Definition: hm_iterative.hh:167
EquationBase::time
TimeGovernor & time()
Definition: equation.hh:148
FieldSet::set_input_list
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
Definition: field_set.hh:289
HM_Iterative::registrar
static const int registrar
Definition: hm_iterative.hh:199
Input::Type::Double
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:534
FLOW123D_FORCE_LINK_IN_CHILD
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
Definition: global_defs.h:104
MPI_DOUBLE
#define MPI_DOUBLE
Definition: mpi.h:156
HM_Iterative::zero_time_step
void zero_time_step() override
Definition: hm_iterative.cc:282
THROW
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Definition: exceptions.hh:53
HM_Iterative::EqFields
Definition: hm_iterative.hh:150
HM_Iterative::EqData
Definition: hm_iterative.hh:136
FieldCommon::flags
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Definition: field_common.hh:191
HM_Iterative::EqFields::ref_pressure_potential
Field< 3, FieldValue< 3 >::Scalar > ref_pressure_potential
Potential of reference (prescribed) pressure from flow b.c. TODO: Swith to BCField when possible.
Definition: hm_iterative.hh:164
fn_fluid_source::operator()
double operator()(double alpha, double beta, double pressure, double old_pressure, double div_u, double old_div_u)
Definition: hm_iterative.cc:63
field_fe.hh
DarcyLMH::get_input_type
static const Input::Type::Record & get_input_type()
Definition: darcy_flow_lmh.cc:104
darcy_flow_lmh.hh
Lumped mixed-hybrid model of linear Darcy flow, possibly unsteady.
FieldCommon::set_time_result_changed
void set_time_result_changed()
Manually mark flag that the field has been changed.
Definition: field_common.hh:734
copy_field
void copy_field(const FieldCommon &from_field_common, FieldFE< dim, Value > &to_field)
Definition: hm_iterative.cc:269
HM_Iterative::EqFields::old_iter_pressure_ptr_
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > old_iter_pressure_ptr_
Definition: hm_iterative.hh:172
fn_hm_coupling_beta
Definition: hm_iterative.cc:42
GenericAssembly::assemble
void assemble(std::shared_ptr< DOFHandlerMultiDim > dh) override
General assemble methods.
Definition: generic_assembly.hh:204
FieldFlag::equation_external_output
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Definition: field_flag.hh:58
HM_Iterative::EqFields::gravity
Field< 3, FieldValue< 3 >::Scalar > gravity
Standard gravity.
Definition: hm_iterative.hh:159
HM_Iterative::eq_data_
EqData eq_data_
Definition: hm_iterative.hh:206
HM_Iterative::update_after_converged
void update_after_converged() override
Save data after iterations have finished.
Definition: hm_iterative.cc:330
FieldCommon::flags_add
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
Definition: field_common.hh:197
HM_Iterative::EqFields::initialize
void initialize(Mesh &mesh, HM_Iterative::EqData &eq_data, const TimeGovernor *time_, double beta_)
Definition: hm_iterative.cc:167
FieldFE::get_dofhandler
std::shared_ptr< DOFHandlerMultiDim > get_dofhandler() const
Definition: field_fe.hh:177
IterativeCoupling::a_tol_
double a_tol_
Absolute tolerance for difference between two succeeding iterations.
Definition: hm_iterative.hh:113
Input::Type::Record::size
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:602
EquationBase::mesh
Mesh & mesh()
Definition: equation.hh:178
Input::Type::Default
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
Input::Type::Record::derive_from
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:196
hm_iterative.hh
fn_fluid_source::time_
const TimeGovernor * time_
Definition: hm_iterative.cc:70
HM_Iterative::HM_Iterative
HM_Iterative(Mesh &mesh, Input::Record in_record)
Definition: hm_iterative.cc:214
HM_Iterative::flow_potential_assembly_
GenericAssembly< FlowPotentialAssemblyHM > * flow_potential_assembly_
Definition: hm_iterative.hh:201
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
FieldFlag::equation_result
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Definition: field_flag.hh:55
IterativeCoupling::record_template
static const Input::Type::Record & record_template()
Definition: hm_iterative.hh:49
sys_profiler.hh
field_model.hh
TimeGovernor
Basic time management functionality for unsteady (and steady) solvers (class Equation).
Definition: time_governor.hh:317
HM_Iterative::update_potential
void update_potential()
Definition: hm_iterative.cc:340
TimeGovernor::step
const TimeStep & step(int index=-1) const
Definition: time_governor.cc:756
fn_fluid_source::fn_fluid_source
fn_fluid_source(const TimeGovernor *time)
Definition: hm_iterative.cc:61
FieldFE::vec
VectorMPI & vec()
Definition: field_fe.hh:181
Input::Type::Default::obligatory
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
HM_Iterative::compute_iteration_error
void compute_iteration_error(double &abs_error, double &rel_error) override
Compute absolute and relative error in the solution.
Definition: hm_iterative.cc:365
FieldCommon
Common abstract parent of all Field<...> classes.
Definition: field_common.hh:76
FieldCommon::print_message_table
static bool print_message_table(ostream &stream, std::string equation_name)
Definition: field_common.cc:96
UnitSI
Class for representation SI units of Fields.
Definition: unit_si.hh:40
assembly_hm.hh
FieldSet::add_coords_field
void add_coords_field()
Definition: field_set.cc:307
EquationBase::input_record_
Input::Record input_record_
Definition: equation.hh:221
Input::Type::Record::declare_key
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
HM_Iterative::EqFields::old_pressure
Field< 3, FieldValue< 3 >::Scalar > old_pressure
Definition: hm_iterative.hh:166
HM_Iterative::update_flow_fields
void update_flow_fields()
Definition: hm_iterative.cc:356
DarcyFlowInterface::get_input_type
static Input::Type::Abstract & get_input_type()
Definition: darcy_flow_interface.hh:23
HM_Iterative::update_solution
void update_solution() override
Definition: hm_iterative.cc:300
TimeGovernor::view
void view(const char *name="") const
Definition: time_governor.cc:772
HM_Iterative::residual_assembly_
GenericAssembly< ResidualAssemblyHM > * residual_assembly_
Definition: hm_iterative.hh:202
Input::Type::Record::close
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:304
HM_Iterative::EqData::p_norm2
double p_norm2
Squared pressure norm in the last iteration.
Definition: hm_iterative.hh:146
Input::Type
Definition: balance.hh:41
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
MPI_SUM
#define MPI_SUM
Definition: mpi.h:196
Field::value
virtual const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Definition: field.hh:447
HM_Iterative::initialize
void initialize() override
Definition: hm_iterative.cc:261
FieldCommon::input_default
FieldCommon & input_default(const string &input_default)
Definition: field_common.hh:139
input_type.hh
FieldFE
Definition: field.hh:59
HM_Iterative::EqData::flow_
std::shared_ptr< DarcyLMH > flow_
steady or unsteady water flow simulator based on MH scheme
Definition: hm_iterative.hh:140
Mesh
Definition: mesh.h:361
HM_Iterative::~HM_Iterative
~HM_Iterative()
Definition: hm_iterative.cc:389
fn_hm_coupling_beta::operator()
double operator()(double alpha, double lame_mu, double lame_lambda, double density, double gravity)
Definition: hm_iterative.cc:47
IterativeCoupling::r_tol_
double r_tol_
Relative tolerance for difference between two succeeding iterations.
Definition: hm_iterative.hh:116
Input::Type::Record::copy_keys
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:216
Input::Type::Array
Class for declaration of inputs sequences.
Definition: type_base.hh:339
Model
Definition: field_model.hh:338
HM_Iterative::EqData::p_dif2
double p_dif2
Squared norm of pressure difference in two subsequent iterations.
Definition: hm_iterative.hh:145
MPI_Allreduce
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
Definition: mpi.h:612
Input::Array
Accessor to input data conforming to declared Array.
Definition: accessors.hh:566
WarningOut
#define WarningOut()
Macro defining 'warning' record of log.
Definition: logger.hh:278
fn_fluid_source
Definition: hm_iterative.cc:59
lame_mu
double lame_mu(double young, double poisson)
Definition: elasticity.cc:79
EquationBase::record_template
static Input::Type::Record & record_template()
Template Record with common keys for derived equations.
Definition: equation.cc:35
MixedPtr
Definition: mixed.hh:247
TimeGovernor::is_default
bool is_default()
Definition: time_governor.hh:388
EquationBase::eq_fieldset_
FieldSet * eq_fieldset_
Definition: equation.hh:228
HM_Iterative::EqData::mechanics_
std::shared_ptr< Elasticity > mechanics_
solute transport with chemistry through operator splitting
Definition: hm_iterative.hh:143
HM_Iterative::update_after_iteration
void update_after_iteration() override
Save data (e.g. solution fields) for the next iteration.
Definition: hm_iterative.cc:322
Elasticity::get_input_type
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
Definition: elasticity.cc:48
HM_Iterative::solve_iteration
void solve_iteration() override
Solve equations and update data (fields).
Definition: hm_iterative.cc:309
lame_lambda
double lame_lambda(double young, double poisson)
Definition: elasticity.cc:85
FieldCommon::description
FieldCommon & description(const string &description)
Definition: field_common.hh:127
Field
Class template representing a field with values dependent on: point, element, and region.
Definition: field.hh:91
fn_hm_coupling_beta::fn_hm_coupling_beta
fn_hm_coupling_beta(double beta_f)
Definition: hm_iterative.cc:44
HM_Iterative::EqFields::old_div_u
Field< 3, FieldValue< 3 >::Scalar > old_div_u
Definition: hm_iterative.hh:168
TimeGovernor::next_time
void next_time()
Proceed to the next time according to current estimated time step.
Definition: time_governor.cc:670
IterativeCoupling::max_it_
unsigned int max_it_
Maximal number of iterations.
Definition: hm_iterative.hh:110
Field::copy_from
void copy_from(const FieldCommon &other) override
Definition: field.impl.hh:367
START_TIMER
#define START_TIMER(tag)
Starts a timer with specified tag.
Definition: sys_profiler.hh:115
GenericAssembly< FlowPotentialAssemblyHM >
DarcyFlowInterface
Definition: darcy_flow_interface.hh:15
FieldCommon::name
FieldCommon & name(const string &name)
Definition: field_common.hh:120
FieldValue
Definition: field_values.hh:645
MessageOut
#define MessageOut()
Macro defining 'message' record of log.
Definition: logger.hh:275