Go to the documentation of this file.
24 #include "petscviewer.h"
25 #include "petscerror.h"
41 *
this += water_content.name(
"water_content")
44 .description(R
"(Water content.
45 It is a fraction of water volume to the whole volume.)");
46 *this += conductivity_richards.name(
"conductivity_richards")
47 .units(
UnitSI().m().s(-1) )
49 .description(
"Computed isotropic scalar conductivity by the soil model.");
51 *
this += water_content_saturated.name(
"water_content_saturated")
52 .description(R
"(Saturated water content (($ \theta_s $)).
53 Relative volume of water in a reference volume of a saturated porous media.)")
57 *
this += water_content_residual.name(
"water_content_residual")
58 .description(R
"(Residual water content (($ \theta_r $)).
59 Relative volume of water in a reference volume of an ideally dry porous media.)")
63 *
this += genuchten_p_head_scale.name(
"genuchten_p_head_scale")
64 .description(R
"(The van Genuchten pressure head scaling parameter (($ \alpha $)).
65 It is related to the inverse of the air entry pressure, i.e. the pressure
66 where the relative water content starts to decrease below 1.)")
70 *
this += genuchten_n_exponent.name(
"genuchten_n_exponent")
71 .description(
"The van Genuchten exponent parameter (($ n $)).")
75 this->set_default_fieldset();
97 auto soil_rec =
it::Record(
"SoilModel",
"Soil model settings.")
100 "Selection of the globally applied soil model. In future we replace this key by a field for selection of the model."
101 "That will allow usage of different soil model in a single simulation.")
103 "Fraction of the water content where we cut and rescale the curve.")
112 "Input data for Darcy flow model.")
114 IT::Default(
"{ \"fields\": [ \"pressure_p0\", \"velocity_p0\" ] }"),
115 "Specification of output fields and output times.")
117 "Soil model settings.")
123 Input::register_class< RichardsLMH, Mesh &, const Input::Record >(
equation_name()) +
147 double fraction= model_rec.val<
double>(
"cut_fraction");
149 eq_data_->soil_model_ = std::make_shared<SoilModel_VanGenuchten>(fraction);
151 eq_data_->soil_model_ = std::make_shared<SoilModel_Irmay>(fraction);
156 eq_data_->water_content_previous_time =
eq_data_->dh_cr_disc_->create_vector();
162 eq_fields_->water_content_ptr = create_field_fe< 3, FieldValue<3>::Scalar >(
eq_data_->dh_cr_disc_);
165 eq_fields_->conductivity_ptr = create_field_fe< 3, FieldValue<3>::Scalar >(
eq_data_->dh_p_);
186 eq_data_->p_edge_solution_previous_time.copy_from(
eq_data_->p_edge_solution);
188 eq_data_->water_content_previous_time.copy_from(water_content_vec);
190 eq_data_->p_edge_solution_previous_time.local_to_ghost_begin();
191 eq_data_->p_edge_solution_previous_time.local_to_ghost_end();
196 return (
eq_fields_->storativity.input_list_size() == 0)
197 && (
eq_fields_->water_content_saturated.input_list_size() == 0);
212 eq_data_->p_edge_solution.local_to_ghost_begin();
213 eq_data_->p_edge_solution.local_to_ghost_end();
232 START_TIMER(
"RichardsLMH::assembly_steady_mh_matrix");
234 END_TIMER(
"RichardsLMH::assembly_steady_mh_matrix");
void accept_time_step() override
postprocess velocity field (add sources)
void read_init_cond_asm() override
Call assemble of read_init_cond_assembly_ and init_cond_postprocess_assembly_.
static UnitSI & dimensionless()
Returns dimensionless unit.
Basic time management class.
std::shared_ptr< FieldSet > eq_fieldset_
const RegionDB & region_db() const
std::shared_ptr< EqData > eq_data_
virtual PetscErrorCode rhs_zero_entries()
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
bool zero_time_term(bool time_global=false) override
virtual PetscErrorCode mat_zero_entries()
static const Input::Type::Record & get_input_type()
GenericAssemblyBase * mh_matrix_assembly_
virtual void assemble(std::shared_ptr< DOFHandlerMultiDim > dh)=0
void assemble(std::shared_ptr< DOFHandlerMultiDim > dh) override
General assemble methods.
void set_matrix_changed()
static const Input::Type::Record & type_field_descriptor()
std::shared_ptr< EqFields > eq_fields_
static const Input::Type::Instance & get_input_type(FieldSet &eq_data, const std::string &equation_name)
void initialize_asm() override
Create and initialize assembly objects.
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Output class for darcy_flow_mh model.
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
static const Input::Type::Record & get_input_type()
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
Basic time management functionality for unsteady (and steady) solvers (class Equation).
Class for representation SI units of Fields.
Input::Record input_record_
RegionSet get_region_set(const std::string &set_name) const
static const std::string field_descriptor_record_description(const string &record_name)
static Input::Type::Abstract & get_input_type()
virtual void start_add_assembly()
std::shared_ptr< EqData > eq_data_
virtual ~RichardsLMH() override
LinSys & lin_sys_schur()
Getter for the linear system of the 2. Schur complement.
virtual void finish_assembly()=0
Global macros to enhance readability and debugging, general constants.
static const int registrar
Registrar of class to factory.
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
GenericAssemblyBase * reconstruct_schur_assembly_
GenericAssembly< InitCondPostprocessAssembly > * init_cond_postprocess_assembly_
general assembly object, hold assembly objects of appropriate dimension
Edge lumped mixed-hybrid solution of unsteady Darcy flow.
RichardsLMH(Mesh &mesh, const Input::Record in_rec, TimeGovernor *tm=nullptr)
GenericAssembly< ReadInitCondAssemblyLMH > * read_init_cond_assembly_
general assembly objects, hold assembly objects of appropriate dimension
void assembly_linear_system() override
#define ASSERT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR) only for debug mode.
std::shared_ptr< EqFields > eq_fields_
#define START_TIMER(tag)
Starts a timer with specified tag.
#define END_TIMER(tag)
Ends a timer with specified tag.
void initialize_specific() override
static std::string equation_name()