Flow123d
JB_transport-112d700
|
Go to the documentation of this file.
35 inline double operator() (
double alpha,
double density,
double gravity,
double pressure)
37 return -alpha*density*gravity*pressure;
63 inline double operator() (
double alpha,
double beta,
double pressure,
double old_pressure,
double div_u,
double old_div_u)
65 return (beta*(pressure-old_pressure) + alpha*(old_div_u - div_u)) /
time_->
dt();
75 std::string equation_name =
"Coupling_Iterative";
77 "Record with data for iterative coupling of flow and mechanics.\n")
84 "Flow equation, provides the velocity field as a result.")
87 "Mechanics, provides the displacement field.")
90 .make_field_descriptor_type(equation_name)),
92 "Input fields of the HM coupling.")
94 "Tuning parameter for iterative splitting. Its default value "
95 "corresponds to a theoretically optimal value with fastest convergence." )
97 "Maximal count of HM iterations." )
99 "Minimal count of HM iterations." )
101 "Absolute tolerance for difference in HM iteration." )
103 "Relative tolerance for difference in HM iteration." )
108 const int HM_Iterative::registrar = Input::register_class< HM_Iterative, Mesh &, const Input::Record >(
"Coupling_Iterative")
121 .
description(
"Volumetric mass density of the fluid.")
127 .
description(
"Gravitational acceleration constant.")
132 *
this +=
beta.
name(
"relaxation_beta")
133 .
description(
"Parameter of numerical method for iterative solution of hydro-mechanical coupling.")
138 .
description(
"Coupling term entering the mechanics equation.")
148 .
description(
"Pressure from last computed iteration.")
153 .
description(
"Displacement divergence from last computed time.")
158 .
description(
"Pressure potential on boundary (taking into account the flow boundary condition.")
163 .
description(
"Coupling term entering the flow equation.")
181 eq_data.
flow_->eq_fields().field_edge_pressure
185 ref_pressure_potential.set(ref_potential_ptr_, 0.0);
196 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);
197 old_pressure.set(old_pressure_ptr_, 0.0);
200 old_iter_pressure.set(old_iter_pressure_ptr_, 0.0);
202 old_div_u_ptr_ = create_field_fe<3, FieldValue<3>::Scalar>(eq_data.
mechanics_->eq_fields().output_div_ptr->get_dofhandler());
203 old_div_u.set(old_div_u_ptr_, 0.0);
209 eq_data.
flow_->eq_fields().field_edge_pressure,
211 eq_data.
mechanics_->eq_fields().output_divergence,
225 using namespace Input;
228 eq_data_ = std::make_shared<EqData>();
244 eq_data_->mechanics_->eq_fields()[
"cross_section"].copy_from(
eq_data_->flow_->eq_fields()[
"cross_section"]);
245 eq_data_->flow_->eq_fields() +=
eq_data_->mechanics_->eq_fields()[
"cross_section_updated"];
246 eq_data_->flow_->eq_fields() +=
eq_data_->mechanics_->eq_fields()[
"stress"];
247 eq_data_->flow_->eq_fields() +=
eq_data_->mechanics_->eq_fields()[
"von_mises_stress"];
248 eq_data_->flow_->eq_fields() +=
eq_data_->mechanics_->eq_fields()[
"mean_stress"];
250 std::stringstream ss;
281 template<
int dim,
class Value>
285 auto vec = to_field.
vec();
289 for (
auto cell : dh->own_range() )
290 vec.set( cell.local_idx(), from_field.
value(cell.elm().centre(), cell.elm()) );
298 std::stringstream ss;
302 eq_data_->mechanics_->update_output_fields();
305 eq_data_->mechanics_->zero_time_step();
309 eq_fields_->old_iter_pressure.set_time_result_changed();
327 eq_data_->flow_->solve_time_step(
false);
331 eq_data_->mechanics_->solve_linear_system();
337 eq_data_->mechanics_->update_output_fields();
339 eq_fields_->old_iter_pressure.set_time_result_changed();
345 eq_data_->flow_->accept_time_step();
347 eq_data_->mechanics_->output_data();
355 auto ref_potential_vec_ =
eq_fields_->ref_potential_ptr_->vec();
356 auto dh =
eq_fields_->ref_potential_ptr_->get_dofhandler();
358 ref_potential_vec_.zero_entries();
361 ref_potential_vec_.local_to_ghost_begin();
362 ref_potential_vec_.local_to_ghost_end();
363 eq_fields_->pressure_potential.set_time_result_changed();
364 eq_fields_->ref_pressure_potential.set_time_result_changed();
372 eq_fields_->flow_source.set_time_result_changed();
380 auto dh =
eq_fields_->old_iter_pressure_ptr_->get_dofhandler();
389 abs_error = sqrt(recv_data[0]);
390 rel_error = abs_error / (sqrt(recv_data[1]) + std::numeric_limits<double>::min());
392 MessageOut().fmt(
"HM Iteration {} abs. difference: {} rel. difference: {}\n"
393 "--------------------------------------------------------",
397 THROW(ExcSolverDiverge() << EI_Reason(
"Reached max_it."));
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
Field< 3, FieldValue< 3 >::Scalar > alpha
Biot coefficient.
static Input::Type::Record & user_fields_template(std::string equation_name)
Template Record with common key user_fields for derived equations.
static const Input::Type::Record & get_input_type()
Define input record.
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Field< 3, FieldValue< 3 >::Scalar > beta
double operator()(double alpha, double density, double gravity, double pressure)
Field< 3, FieldValue< 3 >::Scalar > pressure_potential
Potential -alpha*pressure whose gradient is passed to mechanics as additional load.
const double beta_factor
User-defined factor for iteration parameter.
Field< 3, FieldValue< 3 >::Scalar > density
Density of fluid.
std::shared_ptr< FieldSet > eq_fieldset_
Field< 3, FieldValue< 3 >::Scalar > flow_source
Field< 3, FieldValue< 3 >::Scalar > old_iter_pressure
static const int registrar
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
void zero_time_step() override
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Field< 3, FieldValue< 3 >::Scalar > ref_pressure_potential
Potential of reference (prescribed) pressure from flow b.c. TODO: Swith to BCField when possible.
double operator()(double alpha, double beta, double pressure, double old_pressure, double div_u, double old_div_u)
static const Input::Type::Record & get_input_type()
Lumped mixed-hybrid model of linear Darcy flow, possibly unsteady.
void init_user_fields(Input::Array user_fields, FieldSet &output_fields)
void copy_field(const FieldCommon &from_field_common, FieldFE< dim, Value > &to_field)
void assemble(std::shared_ptr< DOFHandlerMultiDim > dh) override
General assemble methods.
static constexpr Mask equation_external_output
Match an output field, that can be also copy of other field.
Field< 3, FieldValue< 3 >::Scalar > gravity
Standard gravity.
void update_after_converged() override
Save data after iterations have finished.
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
void initialize(Mesh &mesh, HM_Iterative::EqData &eq_data, const TimeGovernor *time_, double beta_)
std::shared_ptr< DOFHandlerMultiDim > get_dofhandler() const
double a_tol_
Absolute tolerance for difference between two succeeding iterations.
const TimeGovernor * time_
void set_default_fieldset()
HM_Iterative(Mesh &mesh, Input::Record in_record)
GenericAssembly< FlowPotentialAssemblyHM > * flow_potential_assembly_
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()
Basic time management functionality for unsteady (and steady) solvers (class Equation).
const TimeStep & step(int index=-1) const
fn_fluid_source(const TimeGovernor *time)
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.
Input::Record input_record_
Field< 3, FieldValue< 3 >::Scalar > old_pressure
Container for various descendants of FieldCommonBase.
void update_flow_fields()
static Input::Type::Abstract & get_input_type()
void update_solution() override
void view(const char *name="") const
GenericAssembly< ResidualAssemblyHM > * residual_assembly_
virtual const Value::return_type & value(const Point &p, const ElementAccessor< spacedim > &elm) const
void initialize() override
FieldCommon & input_default(const string &input_default)
std::shared_ptr< DarcyLMH > flow_
steady or unsteady water flow simulator based on MH scheme
double operator()(double alpha, double lame_mu, double lame_lambda, double density, double gravity)
std::shared_ptr< EqData > eq_data_
double r_tol_
Relative tolerance for difference between two succeeding iterations.
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
#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.
std::shared_ptr< EqFields > eq_fields_
std::shared_ptr< Elasticity > mechanics_
solute transport with chemistry through operator splitting
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.
fn_hm_coupling_beta(double beta_f)
Field< 3, FieldValue< 3 >::Scalar > old_div_u
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.