Flow123d  JS_before_hm-887-g601087d
darcy_flow_lmh.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_lmh.hh
15  * @brief Lumped 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_LMH_HH
34 #define DARCY_FLOW_LMH_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_fe.hh" // for FieldFE
48 #include "fields/field_set.hh" // for FieldSet
49 #include "fields/field_values.hh" // for FieldValue<>::Scalar
50 #include "flow/darcy_flow_interface.hh" // for DarcyFlowInterface
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 #include "la/vector_mpi.hh" // for VectorMPI
59 
60 #include "flow/darcy_flow_mh.hh" // for DarcyMH::EqData
61 
62 class Balance;
63 class DarcyFlowMHOutput;
64 class Element;
65 class Intersection;
66 class LinSys;
67 // class LinSys_BDDC;
68 class LocalSystem;
69 namespace Input {
70  class AbstractRecord;
71  class Record;
72  namespace Type {
73  class Record;
74  class Selection;
75  }
76 }
77 
78 template<int spacedim, class Value> class FieldAddPotential;
79 template<int spacedim, class Value> class FieldDivide;
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 
95 /**
96  * @brief Mixed-hybrid of steady Darcy flow with sources and variable density.
97  *
98  * solve equations:
99  * @f[
100  * q= -{\mathbf{K}} \nabla h -{\mathbf{K}} R \nabla z
101  * @f]
102  * @f[
103  * \mathrm{div} q = f
104  * @f]
105  *
106  * where
107  * - @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.
108  * - @f$ \mathbf{K} @f$ is hydraulic tensor ( its orientation for 2d, 1d case is questionable )
109  * - @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.
110  * Assumes gravity force acts counter to the direction of the @f$ z @f$ axis.
111  * - @f$ R @f$ is destity or gravity variability coefficient. For density driven flow it should be
112  * @f[
113  * R = \frac{\rho}{\rho_0} -1 = \rho_0^{-1}\sum_{i=1}^s c_i
114  * @f]
115  * where @f$ c_i @f$ is concentration in @f$ kg m^{-3} @f$.
116  *
117  * The time key is optional, when not specified the equation is forced to steady regime. Using Steady TimeGovernor which have no dt constraints.
118  *
119  *
120  * TODO:
121  * Make solution regular field (need FeSeystem and parallel DofHandler for edge pressures), then remove get_solution_vector from
122  * Equation interface.
123  */
124 /**
125  * Model for transition coefficients due to Martin, Jaffre, Roberts (see manual for full reference)
126  *
127  * TODO:
128  * - how we can reuse field values computed during assembly
129  *
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 derived from DarcyMH::EqData especially due to the common output class DarcyFlowMHOutput.
147  * This is the only dependence between DarcyMH and DarcyLMH classes.
148  * It is also base class of RichardsLMH::EqData.
149  * */
150  class EqData : public DarcyMH::EqData {
151  public:
152 
153  EqData();
154 
155  std::shared_ptr<SubDOFHandlerMultiDim> dh_p_; ///< DOF handler represents DOFs of element pressure
156 
157  // Propagate test for the time term to the assembly.
158  // This flag is necessary for switching BC to avoid setting zero neumann on the whole boundary in the steady case.
160 
161  // for time term assembly
162  double time_step_;
163 
164  std::shared_ptr<LinSys> lin_sys_schur; //< Linear system of the 2. Schur complement.
165  VectorMPI p_edge_solution; //< 2. Schur complement solution
166  VectorMPI p_edge_solution_previous; //< 2. Schur complement previous solution (iterative)
167  VectorMPI p_edge_solution_previous_time; //< 2. Schur complement previous solution (time)
168 
170  };
171 
172  /// Selection for enum MortarMethod.
173  static const Input::Type::Selection & get_mh_mortar_selection();
174 
175 
176 
177 
178 
179  DarcyLMH(Mesh &mesh, const Input::Record in_rec, TimeGovernor *tm = nullptr);
180 
181  static const Input::Type::Record & type_field_descriptor();
182  static const Input::Type::Record & get_input_type();
183 
184  double last_t() override {
185  return time_->last_t();
186  }
187 
188  std::shared_ptr< FieldFE<3, FieldValue<3>::VectorFixed> > get_velocity_field() override;
189 
190  void init_eq_data();
191  void initialize() override;
192  virtual void initialize_specific();
193  void zero_time_step() override;
194  void update_solution() override;
195 
196  /// Solve the problem without moving to next time and without output.
197  void solve_time_step(bool output = true);
198 
199  /// postprocess velocity field (add sources)
200  virtual void accept_time_step();
201  virtual void postprocess();
202  virtual void output_data() override;
203 
204 
205  EqData &data() { return *data_; }
206 
207  /// Sets external storarivity field (coupling with other equation).
209  { data_->extra_storativity = extra_stor; }
210 
211  /// Sets external source field (coupling with other equation).
212  void set_extra_source(const Field<3, FieldValue<3>::Scalar> &extra_src)
213  { data_->extra_source = extra_src; }
214 
215  virtual ~DarcyLMH() override;
216 
217 
218 protected:
219  /**
220  * Returns true is the fields involved in the time term have values that makes the time term zero.
221  * For time_global==true, it returns true if there are no field descriptors in the input list, so the
222  * fields )of the time ter) have their default values for whole simulation.
223  * If time_global==false (default), only the actual values are considered.
224  */
225  virtual bool zero_time_term(bool time_global=false);
226 
227  /// Solve method common to zero_time_step and update solution.
228  void solve_nonlinear();
229 
230  /**
231  * Create and preallocate MH linear system (including matrix, rhs and solution vectors)
232  */
233  void create_linear_system(Input::AbstractRecord rec);
234 
235  /**
236  * Read initial condition into solution vector.
237  * Must be called after create_linear_system.
238  *
239  * For the LMH scheme we have to be able to save edge pressures in order to
240  * restart simulation or use results of one simulation as initial condition for other one.
241  */
242  void read_initial_condition();
243 
244  /**
245  * In some circumstances, the intial condition must be processed.
246  * It is called at the end of @p read_initial_condition().
247  * This is used in Richards equation due the update of water content.
248  */
249  virtual void initial_condition_postprocess();
250 
251  /**
252  * Allocates linear system matrix for MH.
253  * TODO:
254  * - use general preallocation methods in DofHandler
255  */
256  void allocate_mh_matrix();
257 
258  /**
259  * Assembles linear system matrix for MH.
260  * Element by element assembly is done using dim-template assembly class.
261  * Assembles only steady part of the equation.
262  * TODO:
263  * - include time term - DONE
264  * - add support for Robin type sources
265  * - support for nonlinear solvers - assembly either residual vector, matrix, or both (using FADBAD++)
266  */
267  void assembly_mh_matrix(MultidimAssembly& assembler);
268 
269  void reconstruct_solution_from_schur(MultidimAssembly& assembler);
270 
271  /**
272  * Assembly or update whole linear system.
273  */
274  virtual void assembly_linear_system();
275 
276 // void set_mesh_data_for_bddc(LinSys_BDDC * bddc_ls);
277  /**
278  * Return a norm of residual vector.
279  * TODO: Introduce Equation::compute_residual() updating
280  * residual field, standard part of EqData.
281  */
282  virtual double solution_precision() const;
283 
284  /// Print darcy flow matrix in matlab format into a file.
285  void print_matlab_matrix(string matlab_file);
286 
287  /// Get vector of all DOF indices of given component (0..side, 1..element, 2..edge)
288  std::vector<int> get_component_indices_vec(unsigned int component) const;
289 
290  /// Getter for the linear system of the 2. Schur complement.
292  { return *(data_->lin_sys_schur); }
293 
294  std::shared_ptr<Balance> balance_;
295 
297 
298  int size; // global size of MH matrix
299 
301 
302  // Setting of the nonlinear solver. TODO: Move to the solver class later on.
303  double tolerance_;
304  unsigned int min_n_it_;
305  unsigned int max_n_it_;
306  unsigned int nonlinear_iteration_; //< Actual number of completed nonlinear iterations, need to pass this information into assembly.
307 
308  // Temporary objects holding pointers to appropriate FieldFE
309  // TODO remove after final fix of equations
310  std::shared_ptr<FieldFE<3, FieldValue<3>::VectorFixed>> ele_flux_ptr; ///< Field of flux in barycenter of every element.
311 
312  std::shared_ptr<EqData> data_;
313 
314  friend class DarcyFlowMHOutput;
315  //friend class P0_CouplingAssembler;
316  //friend class P1_CouplingAssembler;
317 
318 private:
319  /// Registrar of class to factory
320  static const int registrar;
321 };
322 
323 #endif //DARCY_FLOW_LMH_HH
324 //-----------------------------------------------------------------------------
325 // vim: set cindent:
326 
VectorMPI p_edge_solution_previous
void set_extra_source(const Field< 3, FieldValue< 3 >::Scalar > &extra_src)
Sets external source field (coupling with other equation).
static const int registrar
Registrar of class to factory.
TYPEDEF_ERR_INFO(EI_KeyName, const string)
VectorMPI p_edge_solution
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
std::shared_ptr< LinSys > lin_sys_schur
Definition: mesh.h:78
std::shared_ptr< SubDOFHandlerMultiDim > dh_p_
DOF handler represents DOFs of element pressure.
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
Basic time management functionality for unsteady (and steady) solvers (class Equation).
Basic time management class.
DarcyFlowMHOutput * output_object
unsigned int min_n_it_
double tolerance_
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
bool data_changed_
unsigned int max_n_it_
std::shared_ptr< EqData > data_
unsigned int nonlinear_iteration_
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:459
void set_extra_storativity(const Field< 3, FieldValue< 3 >::Scalar > &extra_stor)
Sets external storarivity field (coupling with other equation).
EqData & data()
LinSys & lin_sys_schur()
Getter for the linear system of the 2. Schur complement.
std::shared_ptr< Balance > balance_
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)
VectorMPI p_edge_solution_previous_time
mixed-hybrid model of linear Darcy flow, possibly unsteady.
Record type proxy class.
Definition: type_record.hh:182
std::shared_ptr< FieldFE< 3, FieldValue< 3 >::VectorFixed > > ele_flux_ptr
Field of flux in barycenter of every element.
DECLARE_INPUT_EXCEPTION(ExcInputMessage,<< EI_Message::val)
Simple input exception that accepts just string message.
double last_t() override
Return last time of TimeGovernor.
Template for classes storing finite set of named values.
std::map< LongIdx, LocalSystem > seepage_bc_systems