Flow123d  jenkins-Flow123d-linux-release-multijob-282
darcy_flow_mh.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief mixed-hybrid model of linear Darcy flow, possibly unsteady.
27  *
28  * Main object for mixed-hybrid discretization of the linear elliptic PDE (Laplace)
29  * on a multidimensional domain. Discretization of saturated Darcy flow using
30  * RT0 approximation for the velocity
31  *
32  *
33  *
34  * @author Jan Brezina
35  *
36  */
37 
38 /*
39  * list of files dependent on this one:
40  *
41  * posprocess.cc
42  * problem.cc
43  * main.hh
44  * transport.cc
45  */
46 
47 
48 #ifndef DARCY_FLOW_MH_HH
49 #define DARCY_FLOW_MH_HH
50 
51 #include <memory>
52 #include "input/input_type.hh"
53 
54 #include <petscmat.h>
55 #include "system/sys_vector.hh"
56 #include "coupling/equation.hh"
57 #include "flow/mh_dofhandler.hh"
58 #include "input/input_type.hh"
59 #include "la/linsys_BDDC.hh"
60 #include "la/linsys_PETSC.hh"
61 
62 #include "fields/bc_field.hh"
63 #include "fields/field.hh"
64 #include "fields/field_set.hh"
66 #include "flow/old_bcd.hh"
67 
68 
69 /// external types:
70 class LinSys;
71 struct Solver;
72 class Mesh;
73 class SchurComplement;
74 class Distribution;
75 class SparseGraph;
76 class LocalToGlobalMap;
77 class DarcyFlowMHOutput;
78 class Balance;
79 
80 
81 /**
82  * @brief Mixed-hybrid model of linear Darcy flow, possibly unsteady.
83  *
84  * Abstract class for various implementations of Darcy flow. In future there should be
85  * one further level of abstraction for general time dependent problem.
86  *
87  * maybe TODO:
88  * split compute_one_step to :
89  * 1) prepare_next_timestep
90  * 2) actualize_solution - this is for iterative nonlinear solvers
91  *
92  */
93 
94 class DarcyFlowMH : public EquationBase {
95 public:
96  enum MortarMethod {
97  NoMortar = 0,
98  MortarP0 = 1,
100  };
101 
102  /** @brief Data for Darcy flow equation.
103  *
104  */
105  class EqData : public FieldSet {
106  public:
107 
108  /**
109  * For compatibility with old BCD file we have to assign integer codes starting from 1.
110  */
111  enum BC_Type {
112  none=0,
115  robin=3,
117  };
119 
120  /// Collect all fields
121  EqData();
122 
123 
129 
130  BCField<3, FieldValue<3>::Enum > bc_type; // Discrete need Selection for initialization
134 
135  //TODO: these belong to Unsteady flow classes
136  //as long as Unsteady is descendant from Steady, these cannot be transfered..
139 
140  /**
141  * Gravity vector and constant shift of pressure potential. Used to convert piezometric head
142  * to pressure head and vice versa.
143  */
144  arma::vec4 gravity_;
145 
149  };
150 
151 
152  /**
153  * Model for transition coefficients due to Martin, Jaffre, Roberts (see manual for full reference)
154  *
155  * TODO:
156  * - how we can reuse values computed during assembly
157  * we want to make this class see values in
158  *
159  */
161  : EquationBase(mesh, in_rec)
162  {}
163 
169 
170  void get_velocity_seq_vector(Vec &velocity_vec)
171  { velocity_vec = velocity_vector; }
172 
174  double *array;
175  unsigned int size;
176  get_solution_vector(array, size);
177 
178  // here assume that velocity field is extended as constant
179  // to the previous time, so here we set left bound of the interval where the velocity
180  // has current value; this may not be good for every transport !!
181  // we can resolve this when we use FieldFE to store computed velocities in few last steps and
182  // let every equation set time according to nature of the time scheme
183 
184  // in particular this setting is necessary to prevent ConvectinTransport to recreate the transport matrix
185  // every timestep ( this may happen for unsteady flow if we would use time->t() here since it returns infinity.
187  return mh_dh;
188  }
189 
190  virtual void set_concentration_vector(Vec &vc){};
191 
192 
193 protected:
195  double *velocity_array;
196  unsigned int size;
197 
198  get_solution_vector(velocity_array, size);
199  VecCreateSeqWithArray(PETSC_COMM_SELF, 1, mesh_->n_sides(), velocity_array, &velocity_vector);
200 
201  }
202 
203  virtual double solution_precision() const = 0;
204 
207  MH_DofHandler mh_dh; // provides access to seq. solution fluxes and pressures on sides
208 
210 
211  /// object for calculation and writing the water balance to file.
212  boost::shared_ptr<Balance> balance_;
213  /// index of water balance within the Balance object.
214  unsigned int water_balance_idx_;
215 };
216 
217 
218 /**
219  * @brief Mixed-hybrid of steady Darcy flow with sources and variable density.
220  *
221  * solve equations:
222  * @f[
223  * q= -{\mathbf{K}} \nabla h -{\mathbf{K}} R \nabla z
224  * @f]
225  * @f[
226  * \mathrm{div} q = f
227  * @f]
228  *
229  * where
230  * - @f$ q @f$ is flux @f$[ms^{-1}]@f$ for 3d, @f$[m^2s^{-1}]@f$ for 2d and @f$[m^3s^{-1}]@f$ for 1d.
231  * - @f$ \mathbf{K} @f$ is hydraulic tensor ( its orientation for 2d, 1d case is questionable )
232  * - @f$ h = \frac{\pi}{\rho_0 g}+z @f$ is pressure head, @f$ \pi, \rho_0, g @f$ are the pressure, water density, and acceleration of gravity , respectively.
233  * Assumes gravity force acts counter to the direction of the @f$ z @f$ axis.
234  * - @f$ R @f$ is destity or gravity variability coefficient. For density driven flow it should be
235  * @f[
236  * R = \frac{\rho}{\rho_0} -1 = \rho_0^{-1}\sum_{i=1}^s c_i
237  * @f]
238  * where @f$ c_i @f$ is concentration in @f$ kg m^{-3} @f$.
239  *
240  *
241  * TODO:
242  * - consider create auxiliary classes for creation of inverse of block A
243  */
245 {
246 public:
247 
248  class EqData : public DarcyFlowMH::EqData {
249  public:
250 
252  {}
253  };
254 
255  DarcyFlowMH_Steady(Mesh &mesh, const Input::Record in_rec, bool make_tg=true);
256 
258 
259  virtual void update_solution();
260  virtual void get_solution_vector(double * &vec, unsigned int &vec_size);
261  virtual void get_parallel_solution_vector(Vec &vector);
262 
263  /// postprocess velocity field (add sources)
264  virtual void postprocess();
265  virtual void output_data() override;
266 
268 
269 
270 protected:
271  void make_serial_scatter();
272  virtual void modify_system()
273  { ASSERT(0, "Modify system called for Steady darcy.\n"); };
274  virtual void setup_time_term()
275  { ASSERT(0, "Setup time term called for Steady darcy.\n"); };
276 
277 
278  void prepare_parallel( const Input::AbstractRecord in_rec);
279  void make_row_numberings();
280 
281  /**
282  * Create and preallocate MH linear system (including matrix, rhs and solution vectors)
283  */
284  void create_linear_system();
285 
286  /**
287  * Read initial condition into solution vector.
288  * Must be called after create_linear_system.
289  *
290  */
291  virtual void read_init_condition() {};
292 
293  /**
294  * Abstract assembly method used for both assembly and preallocation.
295  * Assembly only steady part of the equation.
296  * TODO:
297  * - use general preallocation methods in DofHandler
298  * - include alos time term
299  * - add support for Robin type sources
300  * - support for nonlinear solvers - assembly either residual vector, matrix, or both (using FADBAD++)
301  *
302  */
304 
305  /// Source term is implemented differently in LMH version.
306  virtual void assembly_source_term();
307 
308  /**
309  * Assembly or update whole linear system.
310  */
311  void assembly_linear_system();
312 
313  void set_mesh_data_for_bddc(LinSys_BDDC * bddc_ls);
314  double solution_precision() const;
315 
316 
318 
319  int size; // global size of MH matrix
320  int n_schur_compls; // number of shur complements to make
321  double *solution; // sequantial scattered solution vector
322 
323 
324  LinSys *schur0; //< whole MH Linear System
325 
326  // parallel
327  Distribution *edge_ds; //< optimal distribution of edges
328  Distribution *el_ds; //< optimal distribution of elements
329  Distribution *side_ds; //< optimal distribution of elements
330  boost::shared_ptr<Distribution> rows_ds; //< final distribution of rows of MH matrix
331 
332  int *el_4_loc; //< array of idexes of local elements (in ordering matching the optimal global)
333  int *row_4_el; //< element index to matrix row
334  int *side_id_4_loc; //< array of ids of local sides
335  int *side_row_4_id; //< side id to matrix row
336  int *edge_4_loc; //< array of indexes of local edges
337  int *row_4_edge; //< edge index to matrix row
338 
339  // MATIS related arrays
340  boost::shared_ptr<LocalToGlobalMap> global_row_4_sub_row; //< global dof index for subdomain index
341 
342  // gather of the solution
343  Vec sol_vec; //< vector over solution array
344  VecScatter par_to_all;
345 
347 
348  friend class DarcyFlowMHOutput;
349  friend class P0_CouplingAssembler;
350  friend class P1_CouplingAssembler;
351 };
352 
353 
355 public:
357  : darcy_(darcy),
358  master_list_(darcy.mesh_->master_elements),
359  intersections_(darcy.mesh_->intersections),
360  master_(nullptr),
361  tensor_average(2)
362  {
363  arma::mat master_map(1,2), slave_map(1,3);
364  master_map.fill(1.0 / 2);
365  slave_map.fill(1.0 / 3);
366 
367  tensor_average[0].push_back( trans( master_map ) * master_map );
368  tensor_average[0].push_back( trans( master_map ) * slave_map );
369  tensor_average[1].push_back( trans( slave_map ) * master_map );
370  tensor_average[1].push_back( trans( slave_map ) * slave_map );
371  }
372 
373  void assembly(LinSys &ls);
374  void pressure_diff(int i_ele,
375  vector<int> &dofs,
376  unsigned int &ele_type,
377  double &delta,
378  arma::vec &dirichlet);
379 private:
381 
383 
386 
388  const Element *master_;
389 
390  /// Row matrices to compute element pressure as average of boundary pressures
392  /// measure of master element, should be sum of intersection measures
393  double delta_0;
394 };
395 
396 
397 
399 public:
401  : darcy_(darcy),
402  intersections_(darcy.mesh_->intersections),
403  rhs(5),
404  dofs(5),
405  dirichlet(5)
406  {
407  rhs.zeros();
408  }
409 
410  void assembly(LinSys &ls);
411  void add_sides(const Element * ele, unsigned int shift, vector<int> &dofs, vector<double> &dirichlet);
412 private:
413 
416 
417  arma::vec rhs;
420 };
421 
422 
423 
424 void mat_count_off_proc_values(Mat m, Vec v);
425 
426 
427 
428 /**
429  * @brief Mixed-hybrid solution of unsteady Darcy flow.
430  *
431  * Standard discretization with time term and sources picewise constant
432  * on the element. This leads to violation of the discrete maximum principle for
433  * non-acute meshes or to too small timesteps. For simplicial meshes this can be solved by lumping to the edges. See DarcyFlowLMH_Unsteady.
434  */
435 
437 {
438 public:
439 
442 
444 protected:
445  void read_init_condition() override;
446  void modify_system() override;
447  void setup_time_term();
448 
449 private:
454 
455 };
456 
457 /**
458  * @brief Edge lumped mixed-hybrid solution of unsteady Darcy flow.
459  *
460  * The time term and sources are evenly distributed form an element to its edges.
461  * This applies directly to the second Schur complement. After this system for pressure traces is solved we reconstruct pressures and side flows as follows:
462  *
463  * -# Element pressure is average of edge pressure. This is in fact same as the MH for steady case so we let SchurComplement class do its job.
464  *
465  * -# We let SchurComplement to reconstruct fluxes and then account time term and sources which are evenly distributed from an element to its sides.
466  * It can be proved, that this keeps continuity of the fluxes over the edges.
467  *
468  * This lumping technique preserves discrete maximum principle for any time step provided one use acute mesh. But in practice even worse meshes are tractable.
469  */
471 {
472 public:
473 
476 
478 protected:
479  void read_init_condition() override;
480  void modify_system() override;
481  void assembly_source_term() override;
482  void setup_time_term();
483  virtual void postprocess();
484 private:
489  //Vec time_term;
490 };
491 
492 #endif //DARCY_FLOW_MH_HH
493 //-----------------------------------------------------------------------------
494 // vim: set cindent:
495 
Abstract base class for equation clasess.
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:51
virtual void postprocess()
postprocess velocity field (add sources)
virtual void get_solution_vector(double *&vec, unsigned int &vec_size)
virtual void setup_time_term()
P1_CouplingAssembler(const DarcyFlowMH_Steady &darcy)
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
vector< vector< arma::mat > > tensor_average
Row matrices to compute element pressure as average of boundary pressures.
MortarMethod mortar_method_
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:52
void mat_count_off_proc_values(Mat m, Vec v)
static Input::Type::Selection mh_mortar_selection
void assembly(LinSys &ls)
bool solution_changed_for_scatter
const DarcyFlowMH_Steady & darcy_
void assembly(LinSys &ls)
Definition: mesh.h:109
static Input::Type::Record input_type
virtual void read_init_condition()
boost::shared_ptr< Distribution > rows_ds
vector< IsecList >::const_iterator ml_it_
BCField< 3, FieldValue< 3 >::Scalar > bc_flux
void pressure_diff(int i_ele, vector< int > &dofs, unsigned int &ele_type, double &delta, arma::vec &dirichlet)
vector< double > dirichlet
const MH_DofHandler & get_mh_dofhandler()
static Input::Type::Record input_type
BCField< 3, FieldValue< 3 >::Enum > bc_type
Field< 3, FieldValue< 3 >::Scalar > storativity
unsigned int n_sides()
Definition: mesh.cc:172
virtual double solution_precision() const =0
class to manage local indices on sub-domain to global indices on domain
void modify_system() override
vector< unsigned int > IsecList
BCField< 3, FieldValue< 3 >::Scalar > bc_pressure
void setup_velocity_vector()
double last_t() const
const DarcyFlowMH_Steady & darcy_
Data for Darcy flow equation.
virtual void update_solution()
Mixed-hybrid solution of unsteady Darcy flow.
virtual void get_parallel_solution_vector(Vec &vector)
virtual void modify_system()
#define ASSERT(...)
Definition: global_defs.h:121
virtual void assembly_source_term()
Source term is implemented differently in LMH version.
boost::shared_ptr< Balance > balance_
object for calculation and writing the water balance to file.
P0_CouplingAssembler(const DarcyFlowMH_Steady &darcy)
Mesh & mesh()
Definition: equation.hh:171
Virtual class for construction and partitioning of a distributed sparse graph.
Definition: sparse_graph.hh:66
Field< 3, FieldValue< 3 >::TensorFixed > anisotropy
void get_velocity_seq_vector(Vec &velocity_vec)
Accessor to the data with type Type::Record.
Definition: accessors.hh:327
void modify_system() override
unsigned int water_balance_idx_
index of water balance within the Balance object.
Class for declaration of polymorphic Record.
Definition: type_record.hh:487
Mesh * mesh_
Definition: equation.hh:218
void read_init_condition() override
const vector< IsecList > & master_list_
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
double solution_precision() const
const vector< Intersection > & intersections_
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:448
const vector< Intersection > & intersections_
virtual void get_solution_vector(double *&vector, unsigned int &size)
Definition: equation.hh:198
Field< 3, FieldValue< 3 >::Scalar > conductivity
static Input::Type::AbstractRecord input_type
double delta_0
measure of master element, should be sum of intersection measures
Distribution * el_ds
Field< 3, FieldValue< 3 >::Scalar > sigma
virtual void postprocess()
postprocess velocity field (add sources)
MH_DofHandler mh_dh
boost::shared_ptr< LocalToGlobalMap > global_row_4_sub_row
virtual void set_concentration_vector(Vec &vc)
void assembly_steady_mh_matrix()
Vector classes to support both Iterator, index and Id access and creating co-located vectors...
Edge lumped mixed-hybrid solution of unsteady Darcy flow.
static std::vector< Input::Type::Record > bc_input_types
Abstract linear system class.
Mixed-hybrid of steady Darcy flow with sources and variable density.
BCField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
static Input::Type::Record input_type
void assembly_source_term() override
Source term is implemented differently in LMH version.
void set_mesh_data_for_bddc(LinSys_BDDC *bddc_ls)
void read_init_condition() override
void add_sides(const Element *ele, unsigned int shift, vector< int > &dofs, vector< double > &dirichlet)
Record type proxy class.
Definition: type_record.hh:169
static Input::Type::Record bc_segment_rec
virtual void output_data() override
Write computed fields.
Field< 3, FieldValue< 3 >::Scalar > cross_section
Distribution * side_ds
Distribution * edge_ds
void set_solution(double time, double *solution, double precision)
DarcyFlowMH(Mesh &mesh, const Input::Record in_rec)
static Input::Type::Selection bc_type_selection
Field< 3, FieldValue< 3 >::Scalar > init_pressure
DarcyFlowMHOutput * output_object
DarcyFlowMH_Steady(Mesh &mesh, const Input::Record in_rec, bool make_tg=true)
CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL.
static Input::Type::AbstractRecord bc_input_type
EqData()
Collect all fields.
Template for classes storing finite set of named values.
Field< 3, FieldValue< 3 >::Scalar > water_source_density
void prepare_parallel(const Input::AbstractRecord in_rec)
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
TimeGovernor * time_
Definition: equation.hh:219
const Element * master_