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")
128 using namespace Input;
134 Record flow_rec = in_record.
val<Record>(
"flow_equation");
138 std::stringstream ss;
143 Record mech_rec = in_record.
val<Record>(
"mechanics_equation");
145 mechanics_->data()[
"cross_section"].copy_from(flow_->data()[
"cross_section"]);
146 mechanics_->initialize();
149 beta_ = in_record.
val<
double>(
"iteration_parameter");
166 template<
int dim,
class Value>
170 auto vec = to_field.
vec();
174 for (
auto cell : dh->own_range() )
175 vec[cell.local_idx()] = from_field.
value(cell.elm().centre(), cell.elm());
183 std::stringstream ss;
187 flow_->zero_time_step();
212 flow_->solve_time_step(
false);
230 flow_->accept_time_step();
231 flow_->output_data();
244 field_edge_pressure.
copy_from(*
flow_->data().field(
"pressure_edge"));
245 for (
auto ele : dh->local_range() )
247 auto elm = ele.elm();
248 LocDofVec dof_indices = ele.get_loc_dof_indices();
249 for (
auto side : ele.side_range() )
254 double pressure = field_edge_pressure.
value(side.centre(), elm);
255 double potential = -alpha*density*gravity*pressure;
257 potential_vec_[dof_indices[side.side_idx()]] = potential;
272 field_ele_pressure.
copy_from(*
flow_->data().field(
"pressure_p0"));
273 for (
auto ele : dh->local_range() )
275 auto elm = ele.elm();
278 double young =
mechanics_->data().young_modulus.value(elm.centre(), elm);
279 double poisson =
mechanics_->data().poisson_ratio.value(elm.centre(), elm);
283 double p = field_ele_pressure.
value(elm.centre(), elm);
286 double src = (beta*(p-old_p) + alpha*(old_div_u - div_u)) /
time_->
dt();
288 beta_vec[ele.local_idx()] = beta;
289 src_vec[ele.local_idx()] = src;
302 double p_dif2 = 0, p_norm2 = 0;
304 field_ele_pressure.
copy_from(*
flow_->data().field(
"pressure_p0"));
305 for (
auto cell : dh->own_range())
307 auto elm = cell.elm();
308 double new_p = field_ele_pressure.
value(elm.centre(), elm);
310 p_dif2 += pow(new_p - old_p, 2)*elm.measure();
311 p_norm2 += pow(old_p, 2)*elm.measure();
314 double send_data[] = { p_dif2, p_norm2 };
317 abs_error = sqrt(recv_data[0]);
318 rel_error = abs_error / sqrt(recv_data[1]);
320 MessageOut().fmt(
"HM Iteration {} abs. difference: {} rel. difference: {}\n" 321 "--------------------------------------------------------",
Common abstract parent of all Field<...> classes.
arma::Col< IntIdx > LocDofVec
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...
double beta_
Tuning parameter for iterative splitting.
static const int registrar
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
#define MessageOut()
Macro defining 'message' record of log.
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...
void initialize(Mesh &mesh)
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.
double a_tol_
Absolute tolerance for difference between two succeeding iterations.
void update_flow_fields()
double lame_mu(double young, double poisson)
void zero_time_step() override
static const Input::Type::Record & record_template()
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > potential_ptr_
FieldFE for pressure_potential field.
const TimeStep & step(int index=-1) const
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.
void update_after_iteration() override
Save data (e.g. solution fields) for the next iteration.
void view(const char *name="") const
static Input::Type::Record & record_template()
Template Record with common keys for derived equations.
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > old_pressure_ptr_
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > div_u_ptr_
void update_solution() override
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
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.
void set_time_result_changed()
Manually mark flag that the field has been changed.
Field< 3, FieldValue< 3 >::Scalar > alpha
Biot coefficient.
std::shared_ptr< Elasticity > mechanics_
solute transport with chemistry through operator splitting
#define START_TIMER(tag)
Starts a timer with specified tag.
virtual Value::return_type const & value(const Point &p, const ElementAccessor< spacedim > &elm) const
Field< 3, FieldValue< 3 >::Scalar > gravity
Standard gravity.
static Input::Type::Abstract & get_input_type()
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
Field< 3, FieldValue< 3 >::Scalar > flow_source
void initialize() override
Field< 3, FieldValue< 3 >::Scalar > beta
std::shared_ptr< DOFHandlerMultiDim > get_dofhandler() const
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
FieldCommon & description(const string &description)
void set_input_list(Input::Array input_list, const TimeGovernor &tg)
static const Input::Type::Record & get_input_type()
void set_field(const RegionSet &domain, FieldBasePtr field, double time=0.0)
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)
void solve_iteration() override
Solve equations and update data (fields).
double lame_lambda(double young, double poisson)
void copy_field(const FieldCommon &from_field_common, FieldFE< dim, Value > &to_field)
#define WarningOut()
Macro defining 'warning' record of log.
FieldCommon & name(const string &name)
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)
FieldCommon & flags(FieldFlag::Flags::Mask mask)
void copy_from(const FieldCommon &other) override
Class for representation SI units of Fields.
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)
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)
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::Scalar > > beta_ptr_