Flow123d  release_3.0.0-859-g84487d0
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 <petscmat.h> // for Mat
37 #include <string.h> // for memcpy
38 #include <algorithm> // for swap
39 #include <boost/exception/info.hpp> // for operator<<, error_inf...
40 #include <memory> // for shared_ptr, allocator...
41 #include <new> // for operator new[]
42 #include <string> // for string, operator<<
43 #include <vector> // for vector, vector<>::con...
44 #include <armadillo>
45 #include "fields/bc_field.hh" // for BCField
46 #include "fields/field.hh" // for Field
47 #include "fields/field_set.hh" // for FieldSet
48 #include "fields/field_values.hh" // for FieldValue<>::Scalar
49 #include "flow/darcy_flow_interface.hh" // for DarcyFlowInterface
50 #include "flow/mh_dofhandler.hh" // for MH_DofHandler, uint
51 #include "input/input_exception.hh" // for DECLARE_INPUT_EXCEPTION
52 #include "input/type_base.hh" // for Array
53 #include "input/type_generic.hh" // for Instance
54 #include "mesh/mesh.h" // for Mesh
55 #include "petscvec.h" // for Vec, _p_Vec, VecScatter
56 #include "system/exceptions.hh" // for ExcStream, operator<<
57 #include "tools/time_governor.hh" // for TimeGovernor
58 
59 class AssemblyBase;
60 class Balance;
61 class DarcyFlowMHOutput;
62 class Element;
63 class Intersection;
64 class LinSys;
65 class LinSys_BDDC;
66 namespace Input {
67  class AbstractRecord;
68  class Record;
69  namespace Type {
70  class Record;
71  class Selection;
72  }
73 }
74 
75 /**
76  * @brief Mixed-hybrid model of linear Darcy flow, possibly unsteady.
77  *
78  * Abstract class for various implementations of Darcy flow. In future there should be
79  * one further level of abstraction for general time dependent problem.
80  *
81  * maybe TODO:
82  * split compute_one_step to :
83  * 1) prepare_next_timestep
84  * 2) actualize_solution - this is for iterative nonlinear solvers
85  *
86  */
87 
88 
89 /**
90  * @brief Mixed-hybrid of steady Darcy flow with sources and variable density.
91  *
92  * solve equations:
93  * @f[
94  * q= -{\mathbf{K}} \nabla h -{\mathbf{K}} R \nabla z
95  * @f]
96  * @f[
97  * \mathrm{div} q = f
98  * @f]
99  *
100  * where
101  * - @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.
102  * - @f$ \mathbf{K} @f$ is hydraulic tensor ( its orientation for 2d, 1d case is questionable )
103  * - @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.
104  * Assumes gravity force acts counter to the direction of the @f$ z @f$ axis.
105  * - @f$ R @f$ is destity or gravity variability coefficient. For density driven flow it should be
106  * @f[
107  * R = \frac{\rho}{\rho_0} -1 = \rho_0^{-1}\sum_{i=1}^s c_i
108  * @f]
109  * where @f$ c_i @f$ is concentration in @f$ kg m^{-3} @f$.
110  *
111  * The time key is optional, when not specified the equation is forced to steady regime. Using Steady TimeGovernor which have no dt constraints.
112  *
113  *
114  * TODO:
115  * Make solution regular field (need FeSeystem and parallel DofHandler for edge pressures), then remove get_solution_vector from
116  * Equation interface.
117  */
118 /**
119  * Model for transition coefficients due to Martin, Jaffre, Roberts (see manual for full reference)
120  *
121  * TODO:
122  * - how we can reuse field values computed during assembly
123  *
124  */
125 
127 {
128 public:
129  TYPEDEF_ERR_INFO( EI_Reason, string);
130  DECLARE_EXCEPTION(ExcSolverDiverge,
131  << "Diverged nonlinear solver. Reason: " << EI_Reason::val
132  );
133  DECLARE_INPUT_EXCEPTION(ExcMissingTimeGovernor,
134  << "Missing the key 'time', obligatory for the transient problems.");
135 
136 
138 
139  /// Type of experimental Mortar-like method for non-compatible 1d-2d interaction.
141  NoMortar = 0,
142  MortarP0 = 1,
143  MortarP1 = 2
144  };
145 
146  /// Class with all fields used in the equation DarcyFlow.
147  /// This is common to all implementations since this provides interface
148  /// to this equation for possible coupling.
149  class EqData : public FieldSet {
150  public:
151 
152  /**
153  * For compatibility with old BCD file we have to assign integer codes starting from 1.
154  */
155  enum BC_Type {
156  none=0,
157  dirichlet=1,
158  total_flux=4,
159  seepage=5,
160  river=6
161  };
162 
163  /// Return a Selection corresponding to enum BC_Type.
164  static const Input::Type::Selection & get_bc_type_selection();
165 
166  /// Creation of all fields.
167  EqData();
168 
169 
175 
176  BCField<3, FieldValue<3>::Enum > bc_type; // Discrete need Selection for initialization
181 
184 
185  /**
186  * Gravity vector and constant shift of pressure potential. Used to convert piezometric head
187  * to pressure head and vice versa.
188  */
189  arma::vec4 gravity_;
191 
192  // Mirroring the following members of DarcyMH:
195 
196 
198 
200 
201  unsigned int local_boundary_index;
202  std::shared_ptr<Balance> balance;
204 
205  unsigned int n_schur_compls;
206  int is_linear; ///< Hack fo BDDC solver.
207  bool force_bc_switch; ///< auxiliary flag for switchting Dirichlet like BC
208 
209  /// Idicator of dirichlet or neumann type of switch boundary conditions.
211  };
212 
213  /// Selection for enum MortarMethod.
214  static const Input::Type::Selection & get_mh_mortar_selection();
215 
216 
217 
218 
219 
220  DarcyMH(Mesh &mesh, const Input::Record in_rec);
221 
222  static const Input::Type::Record & type_field_descriptor();
223  static const Input::Type::Record & get_input_type();
224 
225  const MH_DofHandler &get_mh_dofhandler() override {
226  double *array;
227  unsigned int size;
228  get_solution_vector(array, size);
229 
230  // here assume that velocity field is extended as constant
231  // to the previous time, so here we set left bound of the interval where the velocity
232  // has current value; this may not be good for every transport !!
233  // we can resolve this when we use FieldFE to store computed velocities in few last steps and
234  // let every equation set time according to nature of the time scheme
235 
236  // in particular this setting is necessary to prevent ConvectinTransport to recreate the transport matrix
237  // every timestep ( this may happen for unsteady flow if we would use time->t() here since it returns infinity.
238  mh_dh.set_solution(time_->last_t(), array, solution_precision());
239  return mh_dh;
240  }
241 
242  void init_eq_data();
243  void initialize() override;
244  virtual void initialize_specific();
245  void zero_time_step() override;
246  void update_solution() override;
247 
248  void get_solution_vector(double * &vec, unsigned int &vec_size) override;
249  void get_parallel_solution_vector(Vec &vector) override;
250 
251  /// postprocess velocity field (add sources)
252  virtual void prepare_new_time_step();
253  virtual void postprocess();
254  virtual void output_data() override;
255 
256  virtual ~DarcyMH() override;
257 
258 
259 protected:
260  /**
261  * Returns true is the fields involved in the time term have values that makes the time term zero.
262  * For time_global==true, it returns true if there are no field descriptors in the input list, so the
263  * fields )of the time ter) have their default values for whole simulation.
264  * If time_global==false (default), only the actual values are considered.
265  */
266  virtual bool zero_time_term(bool time_global=false);
267 
268  /// Solve method common to zero_time_step and update solution.
269  void solve_nonlinear();
270  void make_serial_scatter();
271  void modify_system();
272  virtual void setup_time_term();
273 
274 
275  //void prepare_parallel();
276  //void make_row_numberings();
277  /// Initialize global_row_4_sub_row.
278  //void prepare_parallel_bddc();
279 
280  /**
281  * Create and preallocate MH linear system (including matrix, rhs and solution vectors)
282  */
283  void create_linear_system(Input::AbstractRecord rec);
284 
285  /**
286  * Read initial condition into solution vector.
287  * Must be called after create_linear_system.
288  *
289  * For the LMH scheme we have to be able to save edge pressures in order to
290  * restart simulation or use results of one simulation as initial condition for other one.
291  */
292  virtual void read_initial_condition();
293 
294  /**
295  * Part of per element assembly that is specific for MH and LMH respectively.
296  *
297  * This implemnets MH case:
298  * - compute conductivity scaling
299  * - assembly source term
300  * - no time term, managed by diagonal extraction etc.
301  */
302  //virtual void local_assembly_specific(AssemblyData &local_data);
303 
304  /**
305  * Allocates linear system matrix for MH.
306  * TODO:
307  * - use general preallocation methods in DofHandler
308  */
309  void allocate_mh_matrix();
310 
311  /**
312  * Assembles linear system matrix for MH.
313  * Element by element assembly is done using dim-template assembly class.
314  * Assembles only steady part of the equation.
315  * TODO:
316  * - include time term
317  * - add support for Robin type sources
318  * - support for nonlinear solvers - assembly either residual vector, matrix, or both (using FADBAD++)
319  */
320  void assembly_mh_matrix(MultidimAssembly& assembler);
321 
322  /// Source term is implemented differently in LMH version.
323  virtual void assembly_source_term();
324 
325  /**
326  * Assembly or update whole linear system.
327  */
328  virtual void assembly_linear_system();
329 
330  void set_mesh_data_for_bddc(LinSys_BDDC * bddc_ls);
331  /**
332  * Return a norm of residual vector.
333  * TODO: Introduce Equation::compute_residual() updating
334  * residual field, standard part of EqData.
335  */
336  virtual double solution_precision() const;
337 
338  /// Print darcy flow matrix in matlab format into a file.
339  void print_matlab_matrix(string matlab_file);
340 
342  //Vec velocity_vector;
343  MH_DofHandler mh_dh; // provides access to seq. solution fluxes and pressures on sides
344 
345 
346  std::shared_ptr<Balance> balance_;
347 
349 
350  int size; // global size of MH matrix
351  int n_schur_compls; // number of shur complements to make
352  double *solution; // sequantial scattered solution vector
353 
354  // Propagate test for the time term to the assembly.
355  // This flag is necessary for switching BC to avoid setting zero neumann on the whole boundary in the steady case.
358 
359  // Setting of the nonlinear solver. TODO: Move to the solver class later on.
360  double tolerance_;
361  unsigned int min_n_it_;
362  unsigned int max_n_it_;
363  unsigned int nonlinear_iteration_; //< Actual number of completed nonlinear iterations, need to pass this information into assembly.
364 
365 
366  LinSys *schur0; //< whole MH Linear System
367 
368 
369  // gather of the solution
370  Vec sol_vec; //< vector over solution array
371  VecScatter par_to_all;
372 
377 
378  std::shared_ptr<EqData> data_;
379 
380  friend class DarcyFlowMHOutput;
381  //friend class P0_CouplingAssembler;
382  //friend class P1_CouplingAssembler;
383 
384 private:
385  /// Registrar of class to factory
386  static const int registrar;
387 };
388 
389 
390 
391 void mat_count_off_proc_values(Mat m, Vec v);
392 
393 
394 #endif //DARCY_FLOW_MH_HH
395 //-----------------------------------------------------------------------------
396 // vim: set cindent:
397 
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:71
Field< 3, FieldValue< 3 >::Scalar > init_pressure
MH_DofHandler * mh_dh
unsigned int uint
TYPEDEF_ERR_INFO(EI_KeyName, const string)
arma::vec3 gravity_vec_
Abstract linear system class.
Definition: balance.hh:35
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:83
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
Definition: mesh.h:80
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
std::shared_ptr< EqData > data_
void mat_count_off_proc_values(Mat m, Vec v)
Basic time management class.
unsigned int min_n_it_
unsigned int nonlinear_iteration_
BCField< 3, FieldValue< 3 >::Scalar > bc_flux
BCField< 3, FieldValue< 3 >::Scalar > bc_pressure
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_
std::shared_ptr< Balance > balance
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:459
LinSys * schur0
MH_DofHandler mh_dh
MortarMethod
Type of experimental Mortar-like method for non-compatible 1d-2d interaction.
unsigned int local_boundary_index
Vec steady_diagonal
Field< 3, FieldValue< 3 >::TensorFixed > anisotropy
DECLARE_EXCEPTION(ExcWrongDefaultJSON,<< "Consistency Error: Not valid JSON of Default value "<< EI_DefaultStr::qval<< " of type "<< EI_TypeName::qval<< ";\n"<< "During declaration of the key: "<< EI_KeyName::qval)
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
DECLARE_INPUT_EXCEPTION(ExcInputMessage,<< EI_Message::val)
Simple input exception that accepts just string message.
Vec new_diagonal
BCField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
std::vector< std::shared_ptr< AssemblyBase > > MultidimAssembly
Vec previous_solution
Template for classes storing finite set of named values.
MortarMethod mortar_method_
double tolerance_
const MH_DofHandler & get_mh_dofhandler() override
Mixed-hybrid model of linear Darcy flow, possibly unsteady.