Flow123d  release_2.2.0-914-gf1a3a4f
darcy_flow_mh.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file darcy_flow_mh.hh
15  * @brief mixed-hybrid model of linear Darcy flow, possibly unsteady.
16  * @author Jan Brezina
17  *
18  * Main object for mixed-hybrid discretization of the linear elliptic PDE (Laplace)
19  * on a multidimensional domain. Discretization of saturated Darcy flow using
20  * RT0 approximation for the velocity
21  */
22 
23 /*
24  * list of files dependent on this one:
25  *
26  * posprocess.cc
27  * problem.cc
28  * main.hh
29  * transport.cc
30  */
31 
32 
33 #ifndef DARCY_FLOW_MH_HH
34 #define DARCY_FLOW_MH_HH
35 
36 #include <memory>
38 
39 #include <petscmat.h>
40 #include "system/sys_vector.hh"
41 #include "coupling/equation.hh"
42 #include "flow/mh_dofhandler.hh"
43 
44 #include "fields/bc_field.hh"
45 #include "fields/field.hh"
46 #include "fields/field_set.hh"
48 #include "input/input_exception.hh"
49 
50 
51 /// external types:
52 class LinSys;
53 class LinSys_BDDC;
54 class Mesh;
55 class DarcyFlowMHOutput;
56 class Balance;
57 class AssemblyBase;
58 
59 /**
60  * @brief Mixed-hybrid model of linear Darcy flow, possibly unsteady.
61  *
62  * Abstract class for various implementations of Darcy flow. In future there should be
63  * one further level of abstraction for general time dependent problem.
64  *
65  * maybe TODO:
66  * split compute_one_step to :
67  * 1) prepare_next_timestep
68  * 2) actualize_solution - this is for iterative nonlinear solvers
69  *
70  */
71 
72 
73 /**
74  * @brief Mixed-hybrid of steady Darcy flow with sources and variable density.
75  *
76  * solve equations:
77  * @f[
78  * q= -{\mathbf{K}} \nabla h -{\mathbf{K}} R \nabla z
79  * @f]
80  * @f[
81  * \mathrm{div} q = f
82  * @f]
83  *
84  * where
85  * - @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.
86  * - @f$ \mathbf{K} @f$ is hydraulic tensor ( its orientation for 2d, 1d case is questionable )
87  * - @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.
88  * Assumes gravity force acts counter to the direction of the @f$ z @f$ axis.
89  * - @f$ R @f$ is destity or gravity variability coefficient. For density driven flow it should be
90  * @f[
91  * R = \frac{\rho}{\rho_0} -1 = \rho_0^{-1}\sum_{i=1}^s c_i
92  * @f]
93  * where @f$ c_i @f$ is concentration in @f$ kg m^{-3} @f$.
94  *
95  * The time key is optional, when not specified the equation is forced to steady regime. Using Steady TimeGovernor which have no dt constraints.
96  *
97  *
98  * TODO:
99  * Make solution regular field (need FeSeystem and parallel DofHandler for edge pressures), then remove get_solution_vector from
100  * Equation interface.
101  */
102 /**
103  * Model for transition coefficients due to Martin, Jaffre, Roberts (see manual for full reference)
104  *
105  * TODO:
106  * - how we can reuse field values computed during assembly
107  *
108  */
109 
111 {
112 public:
113  TYPEDEF_ERR_INFO( EI_Reason, string);
114  DECLARE_EXCEPTION(ExcSolverDiverge,
115  << "Diverged nonlinear solver. Reason: " << EI_Reason::val
116  );
117  DECLARE_INPUT_EXCEPTION(ExcMissingTimeGovernor,
118  << "Missing the key 'time', obligatory for the transient problems.");
119 
120 
122 
123  /// Type of experimental Mortar-like method for non-compatible 1d-2d interaction.
125  NoMortar = 0,
126  MortarP0 = 1,
128  };
129 
130  /// Class with all fields used in the equation DarcyFlow.
131  /// This is common to all implementations since this provides interface
132  /// to this equation for possible coupling.
133  class EqData : public FieldSet {
134  public:
135 
136  /**
137  * For compatibility with old BCD file we have to assign integer codes starting from 1.
138  */
139  enum BC_Type {
140  none=0,
145  };
146 
147  /// Return a Selection corresponding to enum BC_Type.
149 
150  /// Creation of all fields.
151  EqData();
152 
153 
159 
160  BCField<3, FieldValue<3>::Enum > bc_type; // Discrete need Selection for initialization
165 
168 
169  /**
170  * Gravity vector and constant shift of pressure potential. Used to convert piezometric head
171  * to pressure head and vice versa.
172  */
173  arma::vec4 gravity_;
175 
176  // Mirroring the following members of DarcyMH:
179 
180 
182 
184 
185  unsigned int local_boundary_index;
186  std::shared_ptr<Balance> balance;
188 
189  unsigned int n_schur_compls;
190  int is_linear; ///< Hack fo BDDC solver.
191  bool force_bc_switch; ///< auxiliary flag for switchting Dirichlet like BC
192 
193  /// Idicator of dirichlet or neumann type of switch boundary conditions.
195  };
196 
197  /// Selection for enum MortarMethod.
199 
200 
201 
202 
203 
204  DarcyMH(Mesh &mesh, const Input::Record in_rec);
205 
207  static const Input::Type::Record & get_input_type();
208 
209  const MH_DofHandler &get_mh_dofhandler() override {
210  double *array;
211  unsigned int size;
212  get_solution_vector(array, size);
213 
214  // here assume that velocity field is extended as constant
215  // to the previous time, so here we set left bound of the interval where the velocity
216  // has current value; this may not be good for every transport !!
217  // we can resolve this when we use FieldFE to store computed velocities in few last steps and
218  // let every equation set time according to nature of the time scheme
219 
220  // in particular this setting is necessary to prevent ConvectinTransport to recreate the transport matrix
221  // every timestep ( this may happen for unsteady flow if we would use time->t() here since it returns infinity.
223  return mh_dh;
224  }
225 
226  void init_eq_data();
227  void initialize() override;
228  virtual void initialize_specific();
229  void zero_time_step() override;
230  void update_solution() override;
231 
232  void get_solution_vector(double * &vec, unsigned int &vec_size) override;
233  void get_parallel_solution_vector(Vec &vector) override;
234 
235  /// postprocess velocity field (add sources)
236  virtual void prepare_new_time_step();
237  virtual void postprocess();
238  virtual void output_data() override;
239 
240  virtual ~DarcyMH() override;
241 
242 
243 protected:
244  /**
245  * Returns true is the fields involved in the time term have values that makes the time term zero.
246  * For time_global==true, it returns true if there are no field descriptors in the input list, so the
247  * fields )of the time ter) have their default values for whole simulation.
248  * If time_global==false (default), only the actual values are considered.
249  */
250  virtual bool zero_time_term(bool time_global=false);
251 
252  /// Solve method common to zero_time_step and update solution.
253  void solve_nonlinear();
254  void make_serial_scatter();
255  void modify_system();
256  virtual void setup_time_term();
257 
258 
259  //void prepare_parallel();
260  //void make_row_numberings();
261  /// Initialize global_row_4_sub_row.
262  //void prepare_parallel_bddc();
263 
264  /**
265  * Create and preallocate MH linear system (including matrix, rhs and solution vectors)
266  */
268 
269  /**
270  * Read initial condition into solution vector.
271  * Must be called after create_linear_system.
272  *
273  * For the LMH scheme we have to be able to save edge pressures in order to
274  * restart simulation or use results of one simulation as initial condition for other one.
275  */
276  virtual void read_initial_condition();
277 
278  /**
279  * Part of per element assembly that is specific for MH and LMH respectively.
280  *
281  * This implemnets MH case:
282  * - compute conductivity scaling
283  * - assembly source term
284  * - no time term, managed by diagonal extraction etc.
285  */
286  //virtual void local_assembly_specific(AssemblyData &local_data);
287 
288  /**
289  * Allocates linear system matrix for MH.
290  * TODO:
291  * - use general preallocation methods in DofHandler
292  */
293  void allocate_mh_matrix();
294 
295  /**
296  * Assembles linear system matrix for MH.
297  * Element by element assembly is done using dim-template assembly class.
298  * Assembles only steady part of the equation.
299  * TODO:
300  * - include time term
301  * - add support for Robin type sources
302  * - support for nonlinear solvers - assembly either residual vector, matrix, or both (using FADBAD++)
303  */
304  void assembly_mh_matrix(MultidimAssembly& assembler);
305 
306  /// Source term is implemented differently in LMH version.
307  virtual void assembly_source_term();
308 
309  /**
310  * Assembly or update whole linear system.
311  */
312  virtual void assembly_linear_system();
313 
314  void set_mesh_data_for_bddc(LinSys_BDDC * bddc_ls);
315  /**
316  * Return a norm of residual vector.
317  * TODO: Introduce Equation::compute_residual() updating
318  * residual field, standard part of EqData.
319  */
320  virtual double solution_precision() const;
321 
322  /// Print darcy flow matrix in matlab format into a file.
323  void print_matlab_matrix(string matlab_file);
324 
326  //Vec velocity_vector;
327  MH_DofHandler mh_dh; // provides access to seq. solution fluxes and pressures on sides
328 
329 
330  std::shared_ptr<Balance> balance_;
331 
333 
334  int size; // global size of MH matrix
335  int n_schur_compls; // number of shur complements to make
336  double *solution; // sequantial scattered solution vector
337 
338  // Propagate test for the time term to the assembly.
339  // This flag is necessary for switching BC to avoid setting zero neumann on the whole boundary in the steady case.
342 
343  // Setting of the nonlinear solver. TODO: Move to the solver class later on.
344  double tolerance_;
345  unsigned int min_n_it_;
346  unsigned int max_n_it_;
347  unsigned int nonlinear_iteration_; //< Actual number of completed nonlinear iterations, need to pass this information into assembly.
348 
349 
350  LinSys *schur0; //< whole MH Linear System
351 
352 
353  // gather of the solution
354  Vec sol_vec; //< vector over solution array
355  VecScatter par_to_all;
356 
361 
362  std::shared_ptr<EqData> data_;
363 
364  friend class DarcyFlowMHOutput;
365  //friend class P0_CouplingAssembler;
366  //friend class P1_CouplingAssembler;
367 
368 private:
369  /// Registrar of class to factory
370  static const int registrar;
371 };
372 
373 
374 
375 void mat_count_off_proc_values(Mat m, Vec v);
376 
377 
378 #endif //DARCY_FLOW_MH_HH
379 //-----------------------------------------------------------------------------
380 // vim: set cindent:
381 
void get_solution_vector(double *&vec, unsigned int &vec_size) override
void make_serial_scatter()
Abstract base class for equation clasess.
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:61
Field< 3, FieldValue< 3 >::Scalar > init_pressure
virtual void initialize_specific()
void get_parallel_solution_vector(Vec &vector) override
MH_DofHandler * mh_dh
void assembly_mh_matrix(MultidimAssembly &assembler)
DECLARE_INPUT_EXCEPTION(ExcMissingTimeGovernor,<< "Missing the key 'time', obligatory for the transient problems.")
unsigned int uint
void init_eq_data()
virtual ~DarcyMH() override
virtual void assembly_source_term()
Source term is implemented differently in LMH version.
arma::vec3 gravity_vec_
virtual void output_data() override
Write computed fields.
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:62
static const int registrar
Registrar of class to factory.
unsigned int max_n_it_
VecScatter par_to_all
Field< 3, FieldValue< 3 >::Scalar > sigma
Vec steady_rhs
static const Input::Type::Selection & get_mh_mortar_selection()
Selection for enum MortarMethod.
void initialize() override
Definition: mesh.h:99
static const Input::Type::Record & type_field_descriptor()
arma::vec4 gravity_
std::vector< char > bc_switch_dirichlet
Idicator of dirichlet or neumann type of switch boundary conditions.
int is_linear
Hack fo BDDC solver.
unsigned int n_schur_compls
bool solution_changed_for_scatter
void solve_nonlinear()
Solve method common to zero_time_step and update solution.
std::shared_ptr< EqData > data_
void mat_count_off_proc_values(Mat m, Vec v)
unsigned int min_n_it_
void modify_system()
unsigned int nonlinear_iteration_
double last_t() const
void set_mesh_data_for_bddc(LinSys_BDDC *bddc_ls)
BCField< 3, FieldValue< 3 >::Scalar > bc_flux
void update_solution() override
BCField< 3, FieldValue< 3 >::Scalar > bc_pressure
virtual void postprocess()
Field< 3, FieldValue< 3 >::Scalar > water_source_density
BCField< 3, FieldValue< 3 >::Scalar > bc_switch_pressure
Field< 3, FieldValue< 3 >::Scalar > storativity
std::shared_ptr< Balance > balance_
int n_schur_compls
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
bool force_bc_switch
auxiliary flag for switchting Dirichlet like BC
bool data_changed_
double * solution
bool use_steady_assembly_
void allocate_mh_matrix()
std::shared_ptr< Balance > balance
static const Input::Type::Selection & get_bc_type_selection()
Return a Selection corresponding to enum BC_Type.
virtual void assembly_linear_system()
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:459
LinSys * schur0
DarcyMH(Mesh &mesh, const Input::Record in_rec)
CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL.
void create_linear_system(Input::AbstractRecord rec)
Initialize global_row_4_sub_row.
MH_DofHandler mh_dh
MortarMethod
Type of experimental Mortar-like method for non-compatible 1d-2d interaction.
DECLARE_EXCEPTION(ExcSolverDiverge,<< "Diverged nonlinear solver. Reason: "<< EI_Reason::val)
unsigned int local_boundary_index
Vector classes to support both Iterator, index and Id access and creating co-located vectors...
Vec steady_diagonal
Field< 3, FieldValue< 3 >::TensorFixed > anisotropy
TYPEDEF_ERR_INFO(EI_Reason, string)
EqData()
Creation of all fields.
void print_matlab_matrix(string matlab_file)
Print darcy flow matrix in matlab format into a file.
void zero_time_step() override
virtual void read_initial_condition()
Record type proxy class.
Definition: type_record.hh:182
Field< 3, FieldValue< 3 >::Scalar > cross_section
BCField< 3, FieldValue< 3 >::Enum > bc_type
DarcyFlowMHOutput * output_object
Field< 3, FieldValue< 3 >::Scalar > conductivity
void set_solution(double time, double *solution, double precision)
virtual void prepare_new_time_step()
postprocess velocity field (add sources)
Vec new_diagonal
virtual double solution_precision() const
BCField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
std::vector< std::shared_ptr< AssemblyBase > > MultidimAssembly
virtual void setup_time_term()
Vec previous_solution
Template for classes storing finite set of named values.
MortarMethod mortar_method_
double tolerance_
static const Input::Type::Record & get_input_type()
TimeGovernor * time_
Definition: equation.hh:224
virtual bool zero_time_term(bool time_global=false)
const MH_DofHandler & get_mh_dofhandler() override
Mixed-hybrid model of linear Darcy flow, possibly unsteady.