Flow123d
JS_before_hm-2111-g101b53ca4
|
Go to the documentation of this file.
34 "Record with data for iterative coupling of flow and mechanics.\n")
40 "Flow equation, provides the velocity field as a result.")
42 "Mechanics, provides the displacement field.")
45 .make_field_descriptor_type(
"Coupling_Iterative")),
47 "Input fields of the HM coupling.")
49 "Tuning parameter for iterative splitting. Its default value "
50 "corresponds to a theoretically optimal value with fastest convergence." )
52 "Maximal count of HM iterations." )
54 "Minimal count of HM iterations." )
56 "Absolute tolerance for difference in HM iteration." )
58 "Relative tolerance for difference in HM iteration." )
63 const int HM_Iterative::registrar = Input::register_class< HM_Iterative, Mesh &, const Input::Record >(
"Coupling_Iterative")
76 .
description(
"Volumetric mass density of the fluid.")
82 .
description(
"Gravitational acceleration constant.")
87 *
this +=
beta.
name(
"relaxation_beta")
88 .
description(
"Parameter of numerical method for iterative solution of hydro-mechanical coupling.")
93 .
description(
"Coupling term entering the mechanics equation.")
98 .
description(
"Pressure potential on boundary (taking into account the flow boundary condition.")
103 .
description(
"Coupling term entering the flow equation.")
115 pressure_potential.set(potential_ptr_, 0.0);
118 ref_pressure_potential.set(ref_potential_ptr_, 0.0);
121 beta.set(beta_ptr_, 0.0);
123 flow_source_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(beta_ptr_->get_dofhandler());
124 flow_source.set(flow_source_ptr_, 0.0);
126 old_pressure_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(beta_ptr_->get_dofhandler());
127 old_iter_pressure_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(beta_ptr_->get_dofhandler());
128 div_u_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(beta_ptr_->get_dofhandler());
129 old_div_u_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(beta_ptr_->get_dofhandler());
139 using namespace Input;
149 std::stringstream ss;
156 mechanics_->eq_fields()[
"cross_section"].copy_from(
flow_->eq_fields()[
"cross_section"]);
160 beta_ = in_record.
val<
double>(
"iteration_parameter");
177 template<
int dim,
class Value>
181 auto vec = to_field.
vec();
185 for (
auto cell : dh->own_range() )
186 vec.set( cell.local_idx(), from_field.
value(cell.elm().centre(), cell.elm()) );
194 std::stringstream ss;
198 flow_->zero_time_step();
223 flow_->solve_time_step(
false);
241 flow_->accept_time_step();
242 flow_->output_data();
256 field_edge_pressure.
copy_from(*
flow_->eq_fields().field(
"pressure_edge"));
258 for (
auto ele : dh->local_range() )
260 auto elm = ele.elm();
261 LocDofVec dof_indices = ele.get_loc_dof_indices();
262 for (
auto side : ele.side_range() )
267 double pressure = field_edge_pressure.
value(side.centre(), elm);
268 double potential = -alpha*density*gravity*pressure;
270 potential_vec_.set( dof_indices[side.side_idx()], potential );
274 if (side.side().is_boundary() &&
279 double bc_pressure =
flow_->eq_fields().bc_pressure.value(side.centre(), side.cond().element_accessor());
280 ref_potential_vec_.set(dof_indices[side.side_idx()], -alpha*density*gravity*bc_pressure);
283 ref_potential_vec_.set(dof_indices[side.side_idx()], 0);
287 potential_vec_.local_to_ghost_begin();
288 potential_vec_.local_to_ghost_end();
289 ref_potential_vec_.local_to_ghost_begin();
290 ref_potential_vec_.local_to_ghost_end();
303 field_ele_pressure.
copy_from(*
flow_->eq_fields().field(
"pressure_p0"));
304 for (
auto ele : dh->local_range() )
306 auto elm = ele.elm();
309 double young =
mechanics_->eq_fields().young_modulus.value(elm.centre(), elm);
310 double poisson =
mechanics_->eq_fields().poisson_ratio.value(elm.centre(), elm);
314 double p = field_ele_pressure.
value(elm.centre(), elm);
317 double src = (beta*(p-old_p) + alpha*(old_div_u - div_u)) /
time_->
dt();
319 beta_vec.set(ele.local_idx(), beta);
320 src_vec.set(ele.local_idx(), src);
323 beta_vec.local_to_ghost_begin();
324 src_vec.local_to_ghost_begin();
325 beta_vec.local_to_ghost_end();
326 src_vec.local_to_ghost_end();
337 double p_dif2 = 0, p_norm2 = 0;
339 field_ele_pressure.
copy_from(*
flow_->eq_fields().field(
"pressure_p0"));
340 for (
auto cell : dh->own_range())
342 auto elm = cell.elm();
343 double new_p = field_ele_pressure.
value(elm.centre(), elm);
345 p_dif2 += pow(new_p - old_p, 2)*elm.measure();
346 p_norm2 += pow(old_p, 2)*elm.measure();
349 double send_data[] = { p_dif2, p_norm2 };
352 abs_error = sqrt(recv_data[0]);
353 rel_error = abs_error / (sqrt(recv_data[1]) + std::numeric_limits<double>::min());
355 MessageOut().fmt(
"HM Iteration {} abs. difference: {} rel. difference: {}\n"
356 "--------------------------------------------------------",
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
static const Input::Type::Record & get_input_type()
Define input record.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > old_iter_pressure_ptr_
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
arma::Col< IntIdx > LocDofVec
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > beta_ptr_
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > old_div_u_ptr_
bool set_time(const TimeStep &time, LimitSide limit_side)
std::shared_ptr< RichardsLMH > flow_
steady or unsteady water flow simulator based on MH scheme
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
static const int registrar
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
void zero_time_step() override
Field< 3, FieldValue< 3 >::Scalar > ref_pressure_potential
Potential of reference (prescribed) pressure from flow b.c. TODO: Swith to BCField when possible.
FieldCommon & flags(FieldFlag::Flags::Mask mask)
std::shared_ptr< Elasticity > mechanics_
solute transport with chemistry through operator splitting
void set_time_result_changed()
Manually mark flag that the field has been changed.
void copy_field(const FieldCommon &from_field_common, FieldFE< dim, Value > &to_field)
Field< 3, FieldValue< 3 >::Scalar > pressure_potential
Potential -alpha*pressure whose gradient is passed to mechanics as additional load.
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > flow_source_ptr_
void initialize(Mesh &mesh)
void update_after_converged() override
Save data after iterations have finished.
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > div_u_ptr_
std::shared_ptr< DOFHandlerMultiDim > get_dofhandler() const
double a_tol_
Absolute tolerance for difference between two succeeding iterations.
HM_Iterative(Mesh &mesh, Input::Record in_record)
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > old_pressure_ptr_
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
static const Input::Type::Record & record_template()
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > potential_ptr_
FieldFE for pressure_potential field.
Basic time management functionality for unsteady (and steady) solvers (class Equation).
const TimeStep & step(int index=-1) const
void compute_iteration_error(double &abs_error, double &rel_error) override
Compute absolute and relative error in the solution.
Common abstract parent of all Field<...> classes.
static bool print_message_table(ostream &stream, std::string equation_name)
Class for representation SI units of Fields.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > ref_potential_ptr_
void update_flow_fields()
static Input::Type::Abstract & get_input_type()
void update_solution() override
void view(const char *name="") const
double beta_
Tuning parameter for iterative splitting.
virtual const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm) const
void initialize() override
FieldCommon & input_default(const string &input_default)
Field< 3, FieldValue< 3 >::Scalar > alpha
Biot coefficient.
double r_tol_
Relative tolerance for difference between two succeeding iterations.
static const Input::Type::Record & get_input_type()
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
Field< 3, FieldValue< 3 >::Scalar > gravity
Standard gravity.
#define WarningOut()
Macro defining 'warning' record of log.
double lame_mu(double young, double poisson)
static Input::Type::Record & record_template()
Template Record with common keys for derived equations.
Field< 3, FieldValue< 3 >::Scalar > flow_source
Field< 3, FieldValue< 3 >::Scalar > beta
void update_after_iteration() override
Save data (e.g. solution fields) for the next iteration.
static const Input::Type::Record & get_input_type()
Declare input record type for the equation TransportDG.
void solve_iteration() override
Solve equations and update data (fields).
double lame_lambda(double young, double poisson)
FieldCommon & description(const string &description)
Class template representing a field with values dependent on: point, element, and region.
Field< 3, FieldValue< 3 >::Scalar > density
Density of fluid.
void next_time()
Proceed to the next time according to current estimated time step.
unsigned int max_it_
Maximal number of iterations.
void copy_from(const FieldCommon &other) override
#define START_TIMER(tag)
Starts a timer with specified tag.
FieldCommon & name(const string &name)
#define MessageOut()
Macro defining 'message' record of log.