Flow123d  release_3.0.0-1094-g626f1a1
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 "la/vector_mpi.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  *this += water_content_saturated.name("water_content_saturated")
47  .description(R"(Saturated water content (($ \theta_s $)).
48  Relative volume of water in a reference volume of a saturated porous media.)")
49  .input_default("0.0")
51 
52  *this += water_content_residual.name("water_content_residual")
53  .description(R"(Residual water content (($ \theta_r $)).
54  Relative volume of water in a reference volume of an ideally dry porous media.)")
55  .input_default("0.0")
57 
58  *this += genuchten_p_head_scale.name("genuchten_p_head_scale")
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.)")
62  .input_default("0.0")
63  .units( UnitSI().m(-1) );
64 
65  *this += genuchten_n_exponent.name("genuchten_n_exponent")
66  .description("The van Genuchten exponent parameter (($ n $)).")
67  .input_default("2.0")
69 }
70 
71 
73  it::Record field_descriptor = it::Record("RichardsLMH_Data",FieldCommon::field_descriptor_record_description("RichardsLMH_Data"))
75  .copy_keys( RichardsLMH::EqData().make_field_descriptor_type("RichardsLMH_Data_aux") )
76  .close();
77 
78  auto model_selection = it::Selection("Soil_Model_Type", "")
79  .add_value(SoilModelBase::van_genuchten, "van_genuchten", "Van Genuchten soil model with cutting near zero.")
80  .add_value(SoilModelBase::irmay, "irmay", "Irmay model for conductivity, Van Genuchten model for the water content. Suitable for bentonite.")
81  .close();
82 
83  auto soil_rec = it::Record("SoilModel", "Soil model settings.")
84  .allow_auto_conversion("model_type")
85  .declare_key("model_type", model_selection, it::Default("\"van_genuchten\""),
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.")
88  .declare_key("cut_fraction", it::Double(0.0,1.0), it::Default("0.999"),
89  "Fraction of the water content where we cut and rescale the curve.")
90  .close();
91 
92  return it::Record("Flow_Richards_LMH", "Lumped Mixed-Hybrid solver for unsteady unsaturated Darcy flow.")
95  .declare_key("input_fields", it::Array( field_descriptor ), it::Default::obligatory(),
96  "Input data for Darcy flow model.")
97  .declare_key("soil_model", soil_rec, it::Default("\"van_genuchten\""),
98  "Soil model settings.")
99  .close();
100 }
101 
102 
103 const int RichardsLMH::registrar =
104  Input::register_class< RichardsLMH, Mesh &, const Input::Record >("Flow_Richards_LMH") +
106 
107 
108 
110  : DarcyMH(mesh_in, in_rec)
111 {
112  data_ = make_shared<EqData>();
115  //data_->edge_new_local_4_mesh_idx_ = &(this->edge_new_local_4_mesh_idx_);
116 }
117 
119 
120  auto model_rec = input_record_.val<Input::Record>("soil_model");
121  auto model_type = model_rec.val<SoilModelBase::SoilModelType>("model_type");
122  double fraction= model_rec.val<double>("cut_fraction");
123  if (model_type == SoilModelBase::van_genuchten)
124  data_->soil_model_ = std::make_shared<SoilModel_VanGenuchten>(fraction);
125  else if (model_type == SoilModelBase::irmay)
126  data_->soil_model_ = std::make_shared<SoilModel_Irmay>(fraction);
127  else
128  ASSERT(false);
129 
130  // create edge vectors
131  unsigned int n_local_edges = mh_dh.edge_new_local_4_mesh_idx_.size();
132  unsigned int n_local_sides = mh_dh.side_ds->lsize();
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);
137 
138  Distribution ds_split_edges(n_local_edges, PETSC_COMM_WORLD);
139  vector<int> local_edge_rows(n_local_edges);
140 
141  IS is_loc;
142  for(auto item : mh_dh.edge_new_local_4_mesh_idx_) {
143  local_edge_rows[item.second]=mh_dh.row_4_edge[item.first];
144  }
145  ISCreateGeneral(PETSC_COMM_SELF, local_edge_rows.size(),
146  &(local_edge_rows[0]), PETSC_COPY_VALUES, &(is_loc));
147 
148  VecScatterCreate(schur0->get_solution(), is_loc,
149  data_->phead_edge_.petsc_vec(), PETSC_NULL, &solution_2_edge_scatter_);
150  chkerr(ISDestroy(&is_loc));
151 
152 }
153 
154 
156 {
157 
158 }
159 
160 
162 {
163  // apply initial condition
164  // cycle over local element rows
165  double init_value;
166 
167  for (unsigned int i_loc_el = 0; i_loc_el < mh_dh.el_ds->lsize(); i_loc_el++) {
168  auto ele_ac = mh_dh.accessor(i_loc_el);
169 
170  init_value = data_->init_pressure.value(ele_ac.centre(), ele_ac.element_accessor());
171 
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;
175  VecSetValue(schur0->get_solution(),edge_row, init_value/n_sides_of_edge, ADD_VALUES);
176  }
177  VecSetValue(schur0->get_solution(),ele_ac.ele_row(), init_value,ADD_VALUES);
178  }
179  VecAssemblyBegin(schur0->get_solution());
180  VecAssemblyEnd(schur0->get_solution());
181 
182  // set water_content
183  // pretty ugly since postprocess change fluxes, which cause bad balance, so we must set them back
184  VecCopy(schur0->get_solution(), previous_solution); // store solution vector
185  postprocess();
186  VecSwap(schur0->get_solution(), previous_solution); // restore solution vector
187 
188  //DebugOut() << "init sol:\n";
189  //VecView( schur0->get_solution(), PETSC_VIEWER_STDOUT_WORLD);
190  //DebugOut() << "init water content:\n";
191  //VecView( data_->water_content_previous_it.petsc_vec(), PETSC_VIEWER_STDOUT_WORLD);
192 
194 }
195 
196 
198 {
200  data_->water_content_previous_time.copy(data_->water_content_previous_it);
201  //VecCopy(schur0->get_solution(), previous_solution);
202 }
203 
204 bool RichardsLMH::zero_time_term(bool time_global) {
205  if (time_global) {
206  return (data_->storativity.input_list_size() == 0)
207  && (data_->water_content_saturated.input_list_size() == 0);
208 
209  } else {
210  return (data_->storativity.field_result(mesh_->region_db().get_region_set("BULK"))
211  == result_zeros)
212  && (data_->water_content_saturated.field_result(mesh_->region_db().get_region_set("BULK"))
213  == result_zeros);
214  }
215 }
216 
217 
219 {
220 
221  START_TIMER("RicharsLMH::assembly_linear_system");
222 
223  VecScatterBegin(solution_2_edge_scatter_, schur0->get_solution(), data_->phead_edge_.petsc_vec() , INSERT_VALUES, SCATTER_FORWARD);
224  VecScatterEnd(solution_2_edge_scatter_, schur0->get_solution(), data_->phead_edge_.petsc_vec() , INSERT_VALUES, SCATTER_FORWARD);
225 
226  data_->is_linear = data_->genuchten_p_head_scale.field_result(mesh_->region_db().get_region_set("BULK")) == result_zeros;
227 
228  bool is_steady = zero_time_term();
229  //DebugOut() << "Assembly linear system\n";
230  START_TIMER("full assembly");
231  if (typeid(*schur0) != typeid(LinSys_BDDC)) {
232  schur0->start_add_assembly(); // finish allocation and create matrix
233  }
234  data_->time_step_ = time_->dt();
235  auto multidim_assembler = AssemblyBase::create< AssemblyLMH >(data_);
236 
237 
240 
241  balance_->start_source_assembly(data_->water_balance_idx);
242  balance_->start_mass_assembly(data_->water_balance_idx);
243 
244  assembly_mh_matrix( multidim_assembler ); // fill matrix
245 
246  balance_->finish_source_assembly(data_->water_balance_idx);
247  balance_->finish_mass_assembly(data_->water_balance_idx);
248  //MatView( *const_cast<Mat*>(schur0->get_matrix()), PETSC_VIEWER_STDOUT_WORLD );
249  //VecView( *const_cast<Vec*>(schur0->get_rhs()), PETSC_VIEWER_STDOUT_WORLD);
250 
253 
254 
255  if (! is_steady) {
256  START_TIMER("fix time term");
257  //DebugOut() << "setup time term\n";
258  // assembly time term and rhs
260  }
261 }
262 
263 
264 
266 {
267  FEAL_ASSERT(false).error("Shold not be called.");
268 }
269 
270 
271 
272 
273 
275 
276  // update structures for balance of water volume
278 
279 
280 
281  int side_rows[4];
282  double values[4];
283 
284 
285  VecScatterBegin(solution_2_edge_scatter_, schur0->get_solution(), data_->phead_edge_.petsc_vec() , INSERT_VALUES, SCATTER_FORWARD);
286  VecScatterEnd(solution_2_edge_scatter_, schur0->get_solution(), data_->phead_edge_.petsc_vec() , INSERT_VALUES, SCATTER_FORWARD);
287 
288 
289  // modify side fluxes in parallel
290  // for every local edge take time term on digonal and add it to the corresponding flux
291  //PetscScalar *loc_prev_sol;
292  auto multidim_assembler = AssemblyBase::create< AssemblyLMH >(data_);
293 
294  //VecGetArray(previous_solution, &loc_prev_sol);
295  for (unsigned int i_loc = 0; i_loc < mh_dh.el_ds->lsize(); i_loc++) {
296  auto ele_ac = mh_dh.accessor(i_loc);
297  multidim_assembler[ele_ac.dim()-1]->update_water_content(ele_ac);
298 
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());
302  //double storativity = data_->storativity.value(ele_ac.centre(), ele_ac.element_accessor());
303 
304  for (unsigned int i=0; i<ele_ac.element_accessor()->n_sides(); i++) {
305  //unsigned int loc_edge_row = ele_ac.edge_local_row(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)];
309 
310  values[i] = ele_scale * ele_source - ele_scale * (water_content - water_content_previous_time) / time_->dt();
311  }
312  VecSetValues(schur0->get_solution(), ele_ac.n_sides(), side_rows, values, ADD_VALUES);
313  }
314 
315 
316  VecAssemblyBegin(schur0->get_solution());
317  //VecRestoreArray(previous_solution, &loc_prev_sol);
318  VecAssemblyEnd(schur0->get_solution());
319 
320 }
void initialize_specific() override
FieldSet * eq_data_
Definition: equation.hh:232
Distribution * side_ds
Output class for darcy_flow_mh model.
RegionSet get_region_set(const std::string &set_name) const
Definition: region.cc:329
void assembly_mh_matrix(MultidimAssembly &assembler)
unsigned int uint
unsigned int size() const
Returns number of keys in the Record.
Definition: type_record.hh:598
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:341
virtual PetscErrorCode mat_zero_entries()
Definition: linsys.hh:264
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
Definition: mesh.h:76
static const Input::Type::Record & type_field_descriptor()
Distribution * el_ds
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
Definition: system.hh:147
virtual void finish_assembly()=0
const RegionDB & region_db() const
Definition: mesh.h:143
#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:73
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:346
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:541
void setup_time_term() override
std::shared_ptr< Balance > balance_
FieldCommon & input_default(const string &input_default)
const Vec & get_solution()
Definition: linsys.hh:282
Field< 3, FieldValue< 3 >::Scalar > water_content_saturated
Definition: richards_lmh.hh:71
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:501
static Input::Type::Abstract & get_input_type()
LinSys * schur0
virtual PetscErrorCode rhs_zero_entries()
Definition: linsys.hh:273
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:73
FieldCommon & description(const string &description)
MH_DofHandler mh_dh
LongIdx * row_4_edge
Field< 3, FieldValue< 3 >::Scalar > water_content_residual
Definition: richards_lmh.hh:72
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:72
double dt() const
LocalElementAccessorBase< 3 > accessor(uint local_ele_idx)
Field< 3, FieldValue< 3 >::Scalar > genuchten_n_exponent
Definition: richards_lmh.hh:74
FieldCommon & name(const string &name)
void read_initial_condition() override
void assembly_linear_system() override
void set_matrix_changed()
Definition: linsys.hh:212
Record type proxy class.
Definition: type_record.hh:182
static const int registrar
Registrar of class to factory.
Definition: richards_lmh.hh:99
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:62
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 model of linear Darcy flow, possibly unsteady.
unsigned int lsize(int proc) const
get local size