Flow123d  release_2.2.0-19-g1279b7e
richards_lmh.cc
Go to the documentation of this file.
1 /*
2  * richards_lmh.cc
3  *
4  * Created on: Sep 16, 2015
5  * Author: jb
6  */
7 
8 
9 #include "system/global_defs.h"
10 #include "system/sys_profiler.hh"
11 
12 
13 #include "input/input_type.hh"
14 #include "input/factory.hh"
15 #include "flow/richards_lmh.hh"
17 #include "tools/time_governor.hh"
18 
19 #include "petscmat.h"
20 #include "petscviewer.h"
21 #include "petscerror.h"
22 #include <armadillo>
23 
24 #include "la/schur.hh"
25 
26 #include "coupling/balance.hh"
27 
28 #include "fields/vec_seq_double.hh"
29 
30 // in the third_party/FADBAD++ dir, namespace "fadbad"
31 #include "fadbad.h"
32 #include "badiff.h"
33 #include "fadiff.h"
34 
35 #include "flow/assembly_lmh.hh"
36 
37 
38 FLOW123D_FORCE_LINK_IN_CHILD(richards_lmh);
39 
40 
41 namespace it=Input::Type;
42 
43 
45 {
46  string desc;
47 
48  desc = R"(
49 Saturated water content (($ \theta_s $)).
50 relative volume of the water in a reference volume of a saturated porous media.
51 )" ;
52  ADD_FIELD(water_content_saturated, desc, "0.0");
54 
55  desc = R"(
56 Residual water content (($ \theta_r $)).
57 Relative volume of the water in a reference volume of an ideally dry porous media.
58 )";
59  ADD_FIELD(water_content_residual, desc, "0.0");
61 
62  desc = R"(
63 The van Genuchten pressure head scaling parameter (($ \alpha $)).
64 The parameter of the van Genuchten's model to scale the pressure head.
65 Related to the inverse of the air entry pressure, i.e. the pressure where the relative water content starts to decrease below 1.
66 )";
67  ADD_FIELD(genuchten_p_head_scale, desc, "0.0");
69 
71  "The van Genuchten exponent parameter (($ n $)).", "2.0");
73 
74 }
75 
76 
78  it::Record field_descriptor = it::Record("RichardsLMH_Data",FieldCommon::field_descriptor_record_description("RichardsLMH_Data"))
80  .copy_keys( RichardsLMH::EqData().make_field_descriptor_type("RichardsLMH_Data_aux") )
81  .close();
82 
83  auto model_selection = it::Selection("Soil_Model_Type", "")
84  .add_value(SoilModelBase::van_genuchten, "van_genuchten", "Van Genuchten soil model with cutting near zero.")
85  .add_value(SoilModelBase::irmay, "irmay", "Irmay model for conductivity, Van Genuchten model for the water content. Suitable for bentonite.")
86  .close();
87 
88  auto soil_rec = it::Record("SoilModel", "Setting for the soil model.")
89  .allow_auto_conversion("model_type")
90  .declare_key("model_type", model_selection, it::Default("\"van_genuchten\""),
91  "Selection of the globally applied soil model. In future we replace this key by a field for selection of the model."
92  "That will allow usage of different soil model in a single simulation.")
93  .declare_key("cut_fraction", it::Double(0.0,1.0), it::Default("0.999"),
94  "Fraction of the water content where we cut and rescale the curve.")
95  .close();
96 
97  return it::Record("Flow_Richards_LMH", "Lumped Mixed-Hybrid solver for unsteady saturated Darcy flow.")
100  .declare_key("input_fields", it::Array( field_descriptor ), it::Default::obligatory(),
101  "Input data for Darcy flow model.")
102  .declare_key("soil_model", soil_rec, it::Default("\"van_genuchten\""),
103  "Setting for the soil model.")
104  .close();
105 }
106 
107 
108 const int RichardsLMH::registrar =
109  Input::register_class< RichardsLMH, Mesh &, const Input::Record >("Flow_Richards_LMH") +
111 
112 
113 
115  : DarcyMH(mesh_in, in_rec)
116 {
117  data_ = make_shared<EqData>();
120  //data_->edge_new_local_4_mesh_idx_ = &(this->edge_new_local_4_mesh_idx_);
121 }
122 
124 
125  auto model_rec = input_record_.val<Input::Record>("soil_model");
126  auto model_type = model_rec.val<SoilModelBase::SoilModelType>("model_type");
127  double fraction= model_rec.val<double>("cut_fraction");
128  if (model_type == SoilModelBase::van_genuchten)
129  data_->soil_model_ = std::make_shared<SoilModel_VanGenuchten>(fraction);
130  else if (model_type == SoilModelBase::irmay)
131  data_->soil_model_ = std::make_shared<SoilModel_Irmay>(fraction);
132  else
133  ASSERT(false);
134 
135  // create edge vectors
136  unsigned int n_local_edges = mh_dh.edge_new_local_4_mesh_idx_.size();
137  unsigned int n_local_sides = mh_dh.side_ds->lsize();
138  data_->phead_edge_.resize( n_local_edges);
139  data_->water_content_previous_it.resize(n_local_sides);
140  data_->water_content_previous_time.resize(n_local_sides);
141  data_->capacity.resize(n_local_sides);
142 
143  Distribution ds_split_edges(n_local_edges, PETSC_COMM_WORLD);
144  vector<int> local_edge_rows(n_local_edges);
145 
146  IS is_loc;
147  for(auto item : mh_dh.edge_new_local_4_mesh_idx_) {
148  local_edge_rows[item.second]=mh_dh.row_4_edge[item.first];
149  }
150  ISCreateGeneral(PETSC_COMM_SELF, local_edge_rows.size(),
151  &(local_edge_rows[0]), PETSC_COPY_VALUES, &(is_loc));
152 
153  VecScatterCreate(schur0->get_solution(), is_loc,
154  data_->phead_edge_.petsc_vec(), PETSC_NULL, &solution_2_edge_scatter_);
155  ISDestroy(&is_loc);
156 
157 }
158 
159 
161 {
162 
163 }
164 
165 
167 {
168  // apply initial condition
169  // cycle over local element rows
170  double init_value;
171 
172  for (unsigned int i_loc_el = 0; i_loc_el < mh_dh.el_ds->lsize(); i_loc_el++) {
173  auto ele_ac = mh_dh.accessor(i_loc_el);
174 
175  init_value = data_->init_pressure.value(ele_ac.centre(), ele_ac.element_accessor());
176 
177  FOR_ELEMENT_SIDES(ele_ac.full_iter(),i) {
178  int edge_row = ele_ac.edge_row(i);
179  uint n_sides_of_edge = ele_ac.full_iter()->side(i)->edge()->n_sides;
180  VecSetValue(schur0->get_solution(),edge_row, init_value/n_sides_of_edge, ADD_VALUES);
181  }
182  VecSetValue(schur0->get_solution(),ele_ac.ele_row(), init_value,ADD_VALUES);
183  }
184  VecAssemblyBegin(schur0->get_solution());
185  VecAssemblyEnd(schur0->get_solution());
186 
187  // set water_content
188  // pretty ugly since postprocess change fluxes, which cause bad balance, so we must set them back
189  VecCopy(schur0->get_solution(), previous_solution); // store solution vector
190  postprocess();
191  VecSwap(schur0->get_solution(), previous_solution); // restore solution vector
192 
193  //DebugOut() << "init sol:\n";
194  //VecView( schur0->get_solution(), PETSC_VIEWER_STDOUT_WORLD);
195  //DebugOut() << "init water content:\n";
196  //VecView( data_->water_content_previous_it.petsc_vec(), PETSC_VIEWER_STDOUT_WORLD);
197 
199 }
200 
201 
203 {
205  data_->water_content_previous_time.copy(data_->water_content_previous_it);
206  //VecCopy(schur0->get_solution(), previous_solution);
207 }
208 
209 bool RichardsLMH::zero_time_term(bool time_global) {
210  if (time_global) {
211  return (data_->storativity.input_list_size() == 0)
212  && (data_->water_content_saturated.input_list_size() == 0);
213 
214  } else {
215  return (data_->storativity.field_result(mesh_->region_db().get_region_set("BULK"))
216  == result_zeros)
217  && (data_->water_content_saturated.field_result(mesh_->region_db().get_region_set("BULK"))
218  == result_zeros);
219  }
220 }
221 
222 
224 {
225 
226  START_TIMER("RicharsLMH::assembly_linear_system");
227 
228  VecScatterBegin(solution_2_edge_scatter_, schur0->get_solution(), data_->phead_edge_.petsc_vec() , INSERT_VALUES, SCATTER_FORWARD);
229  VecScatterEnd(solution_2_edge_scatter_, schur0->get_solution(), data_->phead_edge_.petsc_vec() , INSERT_VALUES, SCATTER_FORWARD);
230 
231  is_linear_ = data_->genuchten_p_head_scale.field_result(mesh_->region_db().get_region_set("BULK")) == result_zeros;
232 
233  bool is_steady = zero_time_term();
234  //DebugOut() << "Assembly linear system\n";
235  START_TIMER("full assembly");
236  if (typeid(*schur0) != typeid(LinSys_BDDC)) {
237  schur0->start_add_assembly(); // finish allocation and create matrix
238  }
239  data_->time_step_ = time_->dt();
240  auto multidim_assembler = AssemblyBase::create< AssemblyLMH >(data_);
241 
242 
245 
246  balance_->start_source_assembly(water_balance_idx_);
247  balance_->start_mass_assembly(water_balance_idx_);
248 
249  assembly_mh_matrix( multidim_assembler ); // fill matrix
250 
251  balance_->finish_source_assembly(water_balance_idx_);
252  balance_->finish_mass_assembly(water_balance_idx_);
253  //MatView( *const_cast<Mat*>(schur0->get_matrix()), PETSC_VIEWER_STDOUT_WORLD );
254  //VecView( *const_cast<Vec*>(schur0->get_rhs()), PETSC_VIEWER_STDOUT_WORLD);
255 
258 
259 
260  if (! is_steady) {
261  START_TIMER("fix time term");
262  //DebugOut() << "setup time term\n";
263  // assembly time term and rhs
265  }
266 }
267 
268 
269 
271 {
272  FEAL_ASSERT(false).error("Shold not be called.");
273 }
274 
275 
276 
277 
278 
280 
281  // update structures for balance of water volume
283 
284 
285 
286  int side_rows[4];
287  double values[4];
288 
289 
290  VecScatterBegin(solution_2_edge_scatter_, schur0->get_solution(), data_->phead_edge_.petsc_vec() , INSERT_VALUES, SCATTER_FORWARD);
291  VecScatterEnd(solution_2_edge_scatter_, schur0->get_solution(), data_->phead_edge_.petsc_vec() , INSERT_VALUES, SCATTER_FORWARD);
292 
293 
294  // modify side fluxes in parallel
295  // for every local edge take time term on digonal and add it to the corresponding flux
296  //PetscScalar *loc_prev_sol;
297  auto multidim_assembler = AssemblyBase::create< AssemblyLMH >(data_);
298 
299  //VecGetArray(previous_solution, &loc_prev_sol);
300  for (unsigned int i_loc = 0; i_loc < mh_dh.el_ds->lsize(); i_loc++) {
301  auto ele_ac = mh_dh.accessor(i_loc);
302  multidim_assembler[ele_ac.dim()-1]->update_water_content(ele_ac);
303 
304  double ele_scale = ele_ac.measure() *
305  data_->cross_section.value(ele_ac.centre(), ele_ac.element_accessor()) / ele_ac.n_sides();
306  double ele_source = data_->water_source_density.value(ele_ac.centre(), ele_ac.element_accessor());
307  //double storativity = data_->storativity.value(ele_ac.centre(), ele_ac.element_accessor());
308 
309  FOR_ELEMENT_SIDES(ele_ac.full_iter(),i) {
310  //unsigned int loc_edge_row = ele_ac.edge_local_row(i);
311  side_rows[i] = ele_ac.side_row(i);
312  double water_content = data_->water_content_previous_it[ele_ac.side_local_idx(i)];
313  double water_content_previous_time = data_->water_content_previous_time[ele_ac.side_local_idx(i)];
314 
315  values[i] = ele_scale * ele_source - ele_scale * (water_content - water_content_previous_time) / time_->dt();
316  }
317  VecSetValues(schur0->get_solution(), ele_ac.n_sides(), side_rows, values, ADD_VALUES);
318  }
319 
320 
321  VecAssemblyBegin(schur0->get_solution());
322  //VecRestoreArray(previous_solution, &loc_prev_sol);
323  VecAssemblyEnd(schur0->get_solution());
324 
325 }
void initialize_specific() override
FieldSet * eq_data_
Definition: equation.hh:232
Distribution * side_ds
Output class for darcy_flow_mh model.
void assembly_mh_matrix(MultidimAssembler ma)
unsigned int uint
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:598
int is_linear_
void postprocess() override
Class Input::Type::Default specifies default value of keys of a Input::Type::Record.
Definition: type_record.hh:61
virtual void start_add_assembly()
Definition: linsys.hh:297
virtual PetscErrorCode mat_zero_entries()
Definition: linsys.hh:239
bool zero_time_term(bool time_global=false) override
static Default obligatory()
The factory function to make an empty default value which is obligatory.
Definition: type_record.hh:110
RegionSet get_region_set(const string &set_name) const
Definition: region.cc:328
Definition: mesh.h:97
static const Input::Type::Record & type_field_descriptor()
Distribution * el_ds
#define FOR_ELEMENT_SIDES(i, j)
Definition: elements.h:188
virtual void finish_assembly()=0
const RegionDB & region_db() const
Definition: mesh.h:155
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:346
static const std::string field_descriptor_record_description(const string &record_name)
Definition: field_common.cc:68
#define ADD_FIELD(name,...)
Definition: field_set.hh:269
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
bool solution_changed_for_scatter
std::shared_ptr< EqData > data_
FLOW123D_FORCE_LINK_IN_CHILD(richards_lmh)
Record & close() const
Close the Record for further declarations of keys.
Definition: type_record.cc:303
Basic time management class.
Class for declaration of inputs sequences.
Definition: type_base.hh:345
std::unordered_map< unsigned int, unsigned int > edge_new_local_4_mesh_idx_
void prepare_new_time_step() override
postprocess velocity field (add sources)
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
virtual Record & derive_from(Abstract &parent)
Method to derive new Record from an AbstractRecord parent.
Definition: type_record.cc:195
VecScatter solution_2_edge_scatter_
PETSC scatter from the solution vector to the parallel edge vector with ghost values.
Global macros to enhance readability and debugging, general constants.
virtual Record & allow_auto_conversion(const string &from_key)
Allows shorter input of the Record providing only value of the from_key given as the parameter...
Definition: type_record.cc:132
Class for declaration of the input data that are floating point numbers.
Definition: type_base.hh:540
void setup_time_term() override
std::shared_ptr< Balance > balance_
const Vec & get_solution()
Definition: linsys.hh:257
Field< 3, FieldValue< 3 >::Scalar > water_content_saturated
Definition: richards_lmh.hh:57
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
const Ret val(const string &key) const
#define START_TIMER(tag)
Starts a timer with specified tag.
Mesh * mesh_
Definition: equation.hh:223
Selection & add_value(const int value, const std::string &key, const std::string &description="", TypeBase::attribute_map attributes=TypeBase::attribute_map())
Adds one new value with name given by key to the Selection.
Record & declare_key(const string &key, std::shared_ptr< TypeBase > type, const Default &default_value, const string &description, TypeBase::attribute_map key_attributes=TypeBase::attribute_map())
Declares a new key of the Record.
Definition: type_record.cc:490
static Input::Type::Abstract & get_input_type()
LinSys * schur0
virtual PetscErrorCode rhs_zero_entries()
Definition: linsys.hh:248
std::shared_ptr< EqData > data_
Record & copy_keys(const Record &other)
Copy keys from other record.
Definition: type_record.cc:215
Field< 3, FieldValue< 3 >::Scalar > genuchten_p_head_scale
Definition: richards_lmh.hh:59
MH_DofHandler mh_dh
Field< 3, FieldValue< 3 >::Scalar > water_content_residual
Definition: richards_lmh.hh:58
const Selection & close() const
Close the Selection, no more values can be added.
Input::Record input_record_
Definition: equation.hh:225
static const Input::Type::Record & get_input_type()
Definition: richards_lmh.cc:77
double dt() const
LocalElementAccessorBase< 3 > accessor(uint local_ele_idx)
unsigned int water_balance_idx_
index of water balance within the Balance object.
Field< 3, FieldValue< 3 >::Scalar > genuchten_n_exponent
Definition: richards_lmh.hh:60
void read_initial_condition() override
void assembly_linear_system() override
void set_matrix_changed()
Definition: linsys.hh:187
Record type proxy class.
Definition: type_record.hh:182
static const int registrar
Registrar of class to factory.
Definition: richards_lmh.hh:85
Class for representation SI units of Fields.
Definition: unit_si.hh:40
RichardsLMH(Mesh &mesh, const Input::Record in_rec)
static UnitSI & dimensionless()
Returns dimensionless unit.
Definition: unit_si.cc:55
Vec previous_solution
Template for classes storing finite set of named values.
Input::Type::Record make_field_descriptor_type(const std::string &equation_name) const
Definition: field_set.cc:61
static const Input::Type::Record & get_input_type()
TimeGovernor * time_
Definition: equation.hh:224
void assembly_source_term() override
Source term is implemented differently in LMH version.
#define FEAL_ASSERT(expr)
Definition of assert for debug and release mode.
Definition: asserts.hh:279
Mixed-hybrid of steady Darcy flow with sources and variable density.
unsigned int lsize(int proc) const
get local size