Flow123d  JS_before_hm-915-gc632be9
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 "input/input_exception.hh" // for DECLARE_INPUT_EXCEPTION
51 #include "input/type_base.hh" // for Array
52 #include "input/type_generic.hh" // for Instance
53 #include "mesh/mesh.h" // for Mesh
54 #include "petscvec.h" // for Vec, _p_Vec, VecScatter
55 #include "system/exceptions.hh" // for ExcStream, operator<<
56 #include "tools/time_governor.hh" // for TimeGovernor
57 
58 class AssemblyBase;
59 class Balance;
60 class DarcyFlowMHOutput;
61 class Element;
62 class Intersection;
63 class LinSys;
64 class LinSys_BDDC;
65 namespace Input {
66  class AbstractRecord;
67  class Record;
68  namespace Type {
69  class Record;
70  class Selection;
71  }
72 }
74 
75 template<int spacedim, class Value> class FieldAddPotential;
76 template<int spacedim, class Value> class FieldDivide;
77 
78 /**
79  * @brief Mixed-hybrid model of linear Darcy flow, possibly unsteady.
80  *
81  * Abstract class for various implementations of Darcy flow. In future there should be
82  * one further level of abstraction for general time dependent problem.
83  *
84  * maybe TODO:
85  * split compute_one_step to :
86  * 1) prepare_next_timestep
87  * 2) actualize_solution - this is for iterative nonlinear solvers
88  *
89  */
90 
91 
92 /**
93  * @brief Mixed-hybrid of steady Darcy flow with sources and variable density.
94  *
95  * solve equations:
96  * @f[
97  * q= -{\mathbf{K}} \nabla h -{\mathbf{K}} R \nabla z
98  * @f]
99  * @f[
100  * \mathrm{div} q = f
101  * @f]
102  *
103  * where
104  * - @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.
105  * - @f$ \mathbf{K} @f$ is hydraulic tensor ( its orientation for 2d, 1d case is questionable )
106  * - @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.
107  * Assumes gravity force acts counter to the direction of the @f$ z @f$ axis.
108  * - @f$ R @f$ is destity or gravity variability coefficient. For density driven flow it should be
109  * @f[
110  * R = \frac{\rho}{\rho_0} -1 = \rho_0^{-1}\sum_{i=1}^s c_i
111  * @f]
112  * where @f$ c_i @f$ is concentration in @f$ kg m^{-3} @f$.
113  *
114  * The time key is optional, when not specified the equation is forced to steady regime. Using Steady TimeGovernor which have no dt constraints.
115  *
116  *
117  * TODO:
118  * Make solution regular field (need FeSeystem and parallel DofHandler for edge pressures), then remove get_solution_vector from
119  * Equation interface.
120  */
121 /**
122  * Model for transition coefficients due to Martin, Jaffre, Roberts (see manual for full reference)
123  *
124  * TODO:
125  * - how we can reuse field values computed during assembly
126  *
127  * TODO: Remove in future. It is supposed not to be improved anymore,
128  * however it is kept functioning aside of the LMH lumped version until
129  * the LMH version is stable and optimized.
130  */
131 
133 {
134 public:
135  TYPEDEF_ERR_INFO( EI_Reason, string);
136  DECLARE_EXCEPTION(ExcSolverDiverge,
137  << "Diverged nonlinear solver. Reason: " << EI_Reason::val
138  );
139  DECLARE_INPUT_EXCEPTION(ExcMissingTimeGovernor,
140  << "Missing the key 'time', obligatory for the transient problems.");
141 
142  /** Class with all fields used in the equation DarcyFlow.
143  * This is common to all implementations since this provides interface
144  * to this equation for possible coupling.
145  *
146  * This class is the base class for equation data also in DarcyLMH and RichardsLMH classes
147  * especially due to the common output class DarcyFlowMHOutput.
148  * This is the only dependence between DarcyMH and DarcyLMH classes.
149  */
150  class EqData : public FieldSet {
151  public:
152 
153  /**
154  * For compatibility with old BCD file we have to assign integer codes starting from 1.
155  */
156  enum BC_Type {
157  none=0,
158  dirichlet=1,
159  total_flux=4,
160  seepage=5,
161  river=6
162  };
163 
164  /// Return a Selection corresponding to enum BC_Type.
165  static const Input::Type::Selection & get_bc_type_selection();
166 
167  /// Creation of all fields.
168  EqData();
169 
170 
176 
177  BCField<3, FieldValue<3>::Enum > bc_type; // Discrete need Selection for initialization
182 
185  Field<3, FieldValue<3>::Scalar > extra_storativity; /// Externally added storativity.
186  Field<3, FieldValue<3>::Scalar > extra_source; /// Externally added water source.
187 
192 
193  /**
194  * Gravity vector and constant shift of pressure potential. Used to convert piezometric head
195  * to pressure head and vice versa.
196  */
197  arma::vec4 gravity_;
199 
200  // Mirroring the following members of DarcyMH:
202  std::shared_ptr<DOFHandlerMultiDim> dh_; ///< full DOF handler represents DOFs of sides, elements and edges
203  std::shared_ptr<SubDOFHandlerMultiDim> dh_cr_; ///< DOF handler represents DOFs of edges
204  std::shared_ptr<DOFHandlerMultiDim> dh_cr_disc_; ///< DOF handler represents DOFs of sides
205 
206 
208 
209  MortarMethod mortar_method_;
210 
211  std::shared_ptr<Balance> balance;
213 
214  unsigned int n_schur_compls;
215  int is_linear; ///< Hack fo BDDC solver.
216  bool force_no_neumann_bc; ///< auxiliary flag for switchting Dirichlet like BC
217 
218  /// Idicator of dirichlet or neumann type of switch boundary conditions.
220 
221  VectorMPI full_solution; //< full solution [vel,press,lambda] from 2. Schur complement
222 
224  };
225 
226  /// Selection for enum MortarMethod.
227  static const Input::Type::Selection & get_mh_mortar_selection();
228 
229 
230 
231 
232 
233  DarcyMH(Mesh &mesh, const Input::Record in_rec, TimeGovernor *tm = nullptr);
234 
235  static const Input::Type::Record & type_field_descriptor();
236  static const Input::Type::Record & get_input_type();
237 
238  double last_t() override {
239  return time_->last_t();
240  }
241 
242  std::shared_ptr< FieldFE<3, FieldValue<3>::VectorFixed> > get_velocity_field() override;
243 
244  void init_eq_data();
245  void initialize() override;
246  virtual void initialize_specific();
247  void zero_time_step() override;
248  void update_solution() override;
249  /// Solve the problem without moving to next time and without output.
250  void solve_time_step(bool output = true);
251 
252  /// postprocess velocity field (add sources)
253  virtual void postprocess();
254  virtual void output_data() override;
255 
256  EqData &data() { return *data_; }
257 
259  { data_->extra_storativity = extra_stor; }
260 
261  void set_extra_source(const Field<3, FieldValue<3>::Scalar> &extra_src)
262  { data_->extra_source = extra_src; }
263 
264  virtual ~DarcyMH() override;
265 
266 
267 protected:
268  /**
269  * Returns true is the fields involved in the time term have values that makes the time term zero.
270  * For time_global==true, it returns true if there are no field descriptors in the input list, so the
271  * fields )of the time ter) have their default values for whole simulation.
272  * If time_global==false (default), only the actual values are considered.
273  */
274  virtual bool zero_time_term(bool time_global=false);
275 
276  /// Solve method common to zero_time_step and update solution.
277  void solve_nonlinear();
278  void modify_system();
279  virtual void setup_time_term();
280  void prepare_new_time_step();
281 
282 
283  //void prepare_parallel();
284  //void make_row_numberings();
285 
286  /**
287  * Create and preallocate MH linear system (including matrix, rhs and solution vectors)
288  */
289  void create_linear_system(Input::AbstractRecord rec);
290 
291  /**
292  * Read initial condition into solution vector.
293  * Must be called after create_linear_system.
294  *
295  * For the LMH scheme we have to be able to save edge pressures in order to
296  * restart simulation or use results of one simulation as initial condition for other one.
297  */
298  virtual void read_initial_condition();
299 
300  /**
301  * Part of per element assembly that is specific for MH and LMH respectively.
302  *
303  * This implemnets MH case:
304  * - compute conductivity scaling
305  * - assembly source term
306  * - no time term, managed by diagonal extraction etc.
307  */
308  //virtual void local_assembly_specific(AssemblyData &local_data);
309 
310  /**
311  * Allocates linear system matrix for MH.
312  * TODO:
313  * - use general preallocation methods in DofHandler
314  */
315  void allocate_mh_matrix();
316 
317  /**
318  * Assembles linear system matrix for MH.
319  * Element by element assembly is done using dim-template assembly class.
320  * Assembles only steady part of the equation.
321  * TODO:
322  * - include time term
323  * - add support for Robin type sources
324  * - support for nonlinear solvers - assembly either residual vector, matrix, or both (using FADBAD++)
325  */
326  void assembly_mh_matrix(MultidimAssembly& assembler);
327 
328  /// Source term is implemented differently in LMH version.
329  virtual void assembly_source_term();
330 
331  /**
332  * Assembly or update whole linear system.
333  */
334  virtual void assembly_linear_system();
335 
336  void set_mesh_data_for_bddc(LinSys_BDDC * bddc_ls);
337  /**
338  * Return a norm of residual vector.
339  * TODO: Introduce Equation::compute_residual() updating
340  * residual field, standard part of EqData.
341  */
342  virtual double solution_precision() const;
343 
344  /// Print darcy flow matrix in matlab format into a file.
345  void print_matlab_matrix(string matlab_file);
346 
347  /// Get vector of all DOF indices of given component (0..side, 1..element, 2..edge)
348  std::vector<int> get_component_indices_vec(unsigned int component) const;
349 
350  std::shared_ptr<Balance> balance_;
351 
353 
354  int size; // global size of MH matrix
355  int n_schur_compls; // number of shur complements to make
356 
357  // Propagate test for the time term to the assembly.
358  // This flag is necessary for switching BC to avoid setting zero neumann on the whole boundary in the steady case.
361 
362  // Setting of the nonlinear solver. TODO: Move to the solver class later on.
363  double tolerance_;
364  unsigned int min_n_it_;
365  unsigned int max_n_it_;
366  unsigned int nonlinear_iteration_; //< Actual number of completed nonlinear iterations, need to pass this information into assembly.
367 
368 
369  LinSys *schur0; //< whole MH Linear System
370 
375 
376  // Temporary objects holding pointers to appropriate FieldFE
377  // TODO remove after final fix of equations
378  std::shared_ptr<FieldFE<3, FieldValue<3>::VectorFixed>> ele_flux_ptr; ///< Field of flux in barycenter of every element.
379 
380  std::shared_ptr<EqData> data_;
381 
382  friend class DarcyFlowMHOutput;
383  //friend class P0_CouplingAssembler;
384  //friend class P1_CouplingAssembler;
385 
386 private:
387  /// Registrar of class to factory
388  static const int registrar;
389 };
390 
391 
392 
393 void mat_count_off_proc_values(Mat m, Vec v);
394 /// Helper method fills range (min and max) of given component
395 void dofs_range(unsigned int n_dofs, unsigned int &min, unsigned int &max, unsigned int component);
396 
397 
398 #endif //DARCY_FLOW_MH_HH
399 //-----------------------------------------------------------------------------
400 // vim: set cindent:
401 
Field< 3, FieldValue< 3 >::VectorFixed > field_ele_velocity
void set_extra_storativity(const Field< 3, FieldValue< 3 >::Scalar > &extra_stor)
Field< 3, FieldValue< 3 >::Scalar > extra_source
Externally added storativity.
std::vector< std::shared_ptr< AssemblyBase > > MultidimAssembly
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:74
Field< 3, FieldValue< 3 >::Scalar > init_pressure
unsigned int uint
TYPEDEF_ERR_INFO(EI_KeyName, const string)
arma::vec3 gravity_vec_
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > ele_flux_ptr
Field of flux in barycenter of every element.
std::shared_ptr< DOFHandlerMultiDim > dh_
full DOF handler represents DOFs of sides, elements and edges
Abstract linear system class.
Definition: balance.hh:40
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:92
static const int registrar
Registrar of class to factory.
unsigned int max_n_it_
Field< 3, FieldValue< 3 >::Scalar > extra_storativity
Field< 3, FieldValue< 3 >::Scalar > sigma
Vec steady_rhs
Definition: mesh.h:78
arma::vec4 gravity_
std::shared_ptr< SubDOFHandlerMultiDim > dh_cr_
DOF handler represents DOFs of edges.
std::vector< char > bc_switch_dirichlet
Idicator of dirichlet or neumann type of switch boundary conditions.
int is_linear
Hack fo BDDC solver.
Basic time management functionality for unsteady (and steady) solvers (class Equation).
unsigned int n_schur_compls
void set_extra_source(const Field< 3, FieldValue< 3 >::Scalar > &extra_src)
std::shared_ptr< EqData > data_
void mat_count_off_proc_values(Mat m, Vec v)
Basic time management class.
Field< 3, FieldValue< 3 >::Scalar > field_ele_pressure
Externally added water source.
void dofs_range(unsigned int n_dofs, unsigned int &min, unsigned int &max, unsigned int component)
Helper method fills range (min and max) of given component.
MultidimAssembly multidim_assembler
unsigned int min_n_it_
unsigned int nonlinear_iteration_
BCField< 3, FieldValue< 3 >::Scalar > bc_flux
BCField< 3, FieldValue< 3 >::Scalar > bc_pressure
bool force_no_neumann_bc
auxiliary flag for switchting Dirichlet like BC
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 data_changed_
bool use_steady_assembly_
std::shared_ptr< Balance > balance
VectorMPI full_solution
double last_t() override
Return last time of TimeGovernor.
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:459
LinSys * schur0
std::shared_ptr< DOFHandlerMultiDim > dh_cr_disc_
DOF handler represents DOFs of sides.
Field< 3, FieldValue< 3 >::Scalar > field_ele_piezo_head
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)
EqData & data()
Field< 3, FieldValue< 3 >::Scalar > field_edge_pressure
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
Vec previous_solution
Template for classes storing finite set of named values.
MortarMethod mortar_method_
double tolerance_
Mixed-hybrid model of linear Darcy flow, possibly unsteady.