Flow123d
release_3.0.0-968-gc87a28e79
|
Go to the documentation of this file.
20 #include "petscviewer.h"
21 #include "petscerror.h"
47 .
description(R
"(Saturated water content (($ \theta_s $)).
48 Relative volume of water in a reference volume of a saturated porous media.)")
53 .
description(R
"(Residual water content (($ \theta_r $)).
54 Relative volume of water in a reference volume of an ideally dry porous media.)")
59 .
description(R
"(The van Genuchten pressure head scaling parameter (($ \alpha $)).
60 It is related to the inverse of the air entry pressure, i.e. the pressure
61 where the relative water content starts to decrease below 1.)")
66 .
description(
"The van Genuchten exponent parameter (($ n $)).")
83 auto soil_rec =
it::Record(
"SoilModel",
"Soil model settings.")
86 "Selection of the globally applied soil model. In future we replace this key by a field for selection of the model."
87 "That will allow usage of different soil model in a single simulation.")
89 "Fraction of the water content where we cut and rescale the curve.")
92 return it::Record(
"Flow_Richards_LMH",
"Lumped Mixed-Hybrid solver for unsteady unsaturated Darcy flow.")
96 "Input data for Darcy flow model.")
98 "Soil model settings.")
104 Input::register_class< RichardsLMH, Mesh &, const Input::Record >(
"Flow_Richards_LMH") +
112 data_ = make_shared<EqData>();
122 double fraction= model_rec.val<
double>(
"cut_fraction");
124 data_->soil_model_ = std::make_shared<SoilModel_VanGenuchten>(fraction);
126 data_->soil_model_ = std::make_shared<SoilModel_Irmay>(fraction);
133 data_->phead_edge_.resize( n_local_edges);
134 data_->water_content_previous_it.resize(n_local_sides);
135 data_->water_content_previous_time.resize(n_local_sides);
136 data_->capacity.resize(n_local_sides);
138 Distribution ds_split_edges(n_local_edges, PETSC_COMM_WORLD);
145 ISCreateGeneral(PETSC_COMM_SELF, local_edge_rows.size(),
146 &(local_edge_rows[0]), PETSC_COPY_VALUES, &(is_loc));
150 chkerr(ISDestroy(&is_loc));
167 for (
unsigned int i_loc_el = 0; i_loc_el <
mh_dh.
el_ds->
lsize(); i_loc_el++) {
170 init_value =
data_->init_pressure.value(ele_ac.centre(), ele_ac.element_accessor());
172 for (
unsigned int i=0; i<ele_ac.element_accessor()->n_sides(); i++) {
173 int edge_row = ele_ac.edge_row(i);
174 uint n_sides_of_edge = ele_ac.element_accessor().side(i)->edge()->n_sides;
200 data_->water_content_previous_time.copy(
data_->water_content_previous_it);
206 return (
data_->storativity.input_list_size() == 0)
207 && (
data_->water_content_saturated.input_list_size() == 0);
235 auto multidim_assembler = AssemblyBase::create< AssemblyLMH >(
data_);
292 auto multidim_assembler = AssemblyBase::create< AssemblyLMH >(
data_);
297 multidim_assembler[ele_ac.dim()-1]->update_water_content(ele_ac);
299 double ele_scale = ele_ac.
measure() *
300 data_->cross_section.value(ele_ac.centre(), ele_ac.element_accessor()) / ele_ac.n_sides();
301 double ele_source =
data_->water_source_density.value(ele_ac.centre(), ele_ac.element_accessor());
304 for (
unsigned int i=0; i<ele_ac.element_accessor()->n_sides(); i++) {
306 side_rows[i] = ele_ac.side_row(i);
307 double water_content =
data_->water_content_previous_it[ele_ac.side_local_idx(i)];
308 double water_content_previous_time =
data_->water_content_previous_time[ele_ac.side_local_idx(i)];
310 values[i] = ele_scale * ele_source - ele_scale * (water_content - water_content_previous_time) /
time_->
dt();
312 VecSetValues(
schur0->
get_solution(), ele_ac.n_sides(), side_rows, values, ADD_VALUES);
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
std::shared_ptr< Balance > balance_
const Vec & get_solution()
RichardsLMH(Mesh &mesh, const Input::Record in_rec)
static UnitSI & dimensionless()
Returns dimensionless unit.
unsigned int lsize(int proc) const
get local size
Basic time management class.
Field< 3, FieldValue< 3 >::Scalar > water_content_saturated
Field< 3, FieldValue< 3 >::Scalar > water_content_residual
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
static const Input::Type::Record & type_field_descriptor()
virtual PetscErrorCode rhs_zero_entries()
std::unordered_map< unsigned int, unsigned int > edge_new_local_4_mesh_idx_
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
bool zero_time_term(bool time_global=false) override
virtual PetscErrorCode mat_zero_entries()
void postprocess() override
Field< 3, FieldValue< 3 >::Scalar > genuchten_p_head_scale
void setup_time_term() override
void prepare_new_time_step() override
postprocess velocity field (add sources)
void set_matrix_changed()
FLOW123D_FORCE_LINK_IN_CHILD(richards_lmh)
void assembly_mh_matrix(MultidimAssembly &assembler)
#define FEAL_ASSERT(expr)
Definition of assert for debug and release mode.
Output class for darcy_flow_mh model.
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
LocalElementAccessorBase< 3 > accessor(uint local_ele_idx)
std::shared_ptr< EqData > data_
Class for representation SI units of Fields.
const RegionDB & region_db() const
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)
VecScatter solution_2_edge_scatter_
PETSC scatter from the solution vector to the parallel edge vector with ghost values.
static Input::Type::Abstract & get_input_type()
virtual void start_add_assembly()
std::shared_ptr< EqData > data_
FieldCommon & input_default(const string &input_default)
static const Input::Type::Record & get_input_type()
static const Input::Type::Record & get_input_type()
virtual void finish_assembly()=0
Global macros to enhance readability and debugging, general constants.
static const int registrar
Registrar of class to factory.
Field< 3, FieldValue< 3 >::Scalar > genuchten_n_exponent
void assembly_source_term() override
Source term is implemented differently in LMH version.
FieldCommon & description(const string &description)
bool solution_changed_for_scatter
void read_initial_condition() override
void assembly_linear_system() override
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
#define START_TIMER(tag)
Starts a timer with specified tag.
void initialize_specific() override
FieldCommon & name(const string &name)