Flow123d  3.9.1-c3f8cb5
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 
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 #include "la/vector_mpi.hh" // for VectorMPI
58 
59 #include "flow/darcy_flow_mh.hh" // for DarcyMH::EqData
60 
61 class Balance;
62 class DarcyFlowMHOutput;
63 class Element;
64 class Intersection;
65 class LinSys;
66 // class LinSys_BDDC;
67 class LocalSystem;
68 class LocalConstraint;
69 namespace Input {
70  class AbstractRecord;
71  class Record;
72  namespace Type {
73  class Record;
74  class Selection;
75  }
76 }
77 template<unsigned int dim> class ReadInitCondAssemblyLMH;
79 template< template<IntDim...> class DimAssembly> class GenericAssembly;
80 
81 
82 /**
83  * @brief Mixed-hybrid model of linear Darcy flow, possibly unsteady.
84  *
85  * Abstract class for various implementations of Darcy flow. In future there should be
86  * one further level of abstraction for general time dependent problem.
87  *
88  * maybe TODO:
89  * split compute_one_step to :
90  * 1) prepare_next_timestep
91  * 2) actualize_solution - this is for iterative nonlinear solvers
92  *
93  */
94 
95 
96 /**
97  * @brief Mixed-hybrid of steady Darcy flow with sources and variable density.
98  *
99  * solve equations:
100  * @f[
101  * q= -{\mathbf{K}} \nabla h -{\mathbf{K}} R \nabla z
102  * @f]
103  * @f[
104  * \mathrm{div} q = f
105  * @f]
106  *
107  * where
108  * - @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.
109  * - @f$ \mathbf{K} @f$ is hydraulic tensor ( its orientation for 2d, 1d case is questionable )
110  * - @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.
111  * Assumes gravity force acts counter to the direction of the @f$ z @f$ axis.
112  * - @f$ R @f$ is destity or gravity variability coefficient. For density driven flow it should be
113  * @f[
114  * R = \frac{\rho}{\rho_0} -1 = \rho_0^{-1}\sum_{i=1}^s c_i
115  * @f]
116  * where @f$ c_i @f$ is concentration in @f$ kg m^{-3} @f$.
117  *
118  * The time key is optional, when not specified the equation is forced to steady regime. Using Steady TimeGovernor which have no dt constraints.
119  *
120  *
121  * TODO:
122  * Make solution regular field (need FeSeystem and parallel DofHandler for edge pressures), then remove get_solution_vector from
123  * Equation interface.
124  */
125 /**
126  * Model for transition coefficients due to Martin, Jaffre, Roberts (see manual for full reference)
127  *
128  * TODO:
129  * - how we can reuse field values computed during assembly
130  *
131  */
132 
134 {
135 public:
136  TYPEDEF_ERR_INFO( EI_Reason, string);
137  DECLARE_EXCEPTION(ExcSolverDiverge,
138  << "Diverged nonlinear solver. Reason: " << EI_Reason::val
139  );
140  DECLARE_INPUT_EXCEPTION(ExcMissingTimeGovernor,
141  << "Missing the key 'time', obligatory for the transient problems.");
142 
143  /** Class with all fields used in the equation DarcyFlow.
144  * This is common to all implementations since this provides interface
145  * to this equation for possible coupling.
146  *
147  * This class is derived from DarcyMH::EqData especially due to the common output class DarcyFlowMHOutput.
148  * This is the only dependence between DarcyMH and DarcyLMH classes.
149  * It is also base class of RichardsLMH::EqData.
150  * */
151  class EqFields : public DarcyMH::EqFields {
152  public:
153 
154  EqFields();
155 
156  };
157 
158  class EqData : public DarcyMH::EqData {
159  public:
160 
161  EqData();
162 
163  void init(); ///< Initialize vectors, ...
164  void reset(); ///< Reset data members
165 
166  std::shared_ptr<SubDOFHandlerMultiDim> dh_p_; ///< DOF handler represents DOFs of element pressure
167 
168  // Propagate test for the time term to the assembly.
169  // This flag is necessary for switching BC to avoid setting zero neumann on the whole boundary in the steady case.
171 
172  // for time term assembly
173  double time_step_;
174 
175  std::shared_ptr<LinSys> lin_sys_schur; //< Linear system of the 2. Schur complement.
176  VectorMPI p_edge_solution; //< 2. Schur complement solution
177  VectorMPI p_edge_solution_previous; //< 2. Schur complement previous solution (iterative)
178  VectorMPI p_edge_solution_previous_time; //< 2. Schur complement previous solution (time)
179 
181 
182  /// Shared Balance object
183  std::shared_ptr<Balance> balance_;
184 
185  unsigned int nonlinear_iteration_; //< Actual number of completed nonlinear iterations, need to pass this information into assembly.
186 
187  /// Following data members are stored in vectors, one item for every cell.
188 // /** TODO: Investigate why the hell do we need this flag.
189 // * If removed, it does not break any of the integration tests,
190 // * however it must influence the Dirichlet rows in matrix.
191 // */
192 // std::vector<unsigned int> dirichlet_edge;
193 
197  std::array<std::vector<unsigned int>, 3> loc_side_dofs;
198  std::array<std::vector<unsigned int>, 3> loc_edge_dofs;
199  std::array<unsigned int, 3> loc_ele_dof;
200 
201 // // std::shared_ptr<MortarAssemblyBase> mortar_assembly;
202 
203  std::vector<bool> save_local_system_; ///< Flag for saving the local system. Currently used only in case of seepage BC.
204  std::vector<bool> bc_fluxes_reconstruted; ///< Flag indicating whether the fluxes for seepage BC has been reconstructed already.
205  std::array<unsigned int, 3> schur_offset_; ///< Index offset in the local system for the Schur complement (of dim = 1,2,3).
206  };
207 
208  /// Selection for enum MortarMethod.
210 
211 
212 
213 
214 
215  DarcyLMH(Mesh &mesh, const Input::Record in_rec, TimeGovernor *tm = nullptr);
216 
218  static const Input::Type::Record & get_input_type();
219 
220  void init_eq_data();
221  void initialize() override;
222  virtual void initialize_specific();
223  void zero_time_step() override;
224  void update_solution() override;
225 
226  /// Solve the problem without moving to next time and without output.
227  void solve_time_step(bool output = true);
228 
229  /// postprocess velocity field (add sources)
230  virtual void accept_time_step();
231  virtual void postprocess();
232  virtual void output_data() override;
233 
234  virtual double solved_time() override;
235 
236  inline EqFields &eq_fields() { return *eq_fields_; }
237  inline EqData &eq_data() { return *eq_data_; }
238 
239  /// Sets external storarivity field (coupling with other equation).
241  { eq_fields_->extra_storativity = extra_stor; }
242 
243  /// Sets external source field (coupling with other equation).
244  void set_extra_source(const Field<3, FieldValue<3>::Scalar> &extra_src)
245  { eq_fields_->extra_source = extra_src; }
246 
247  virtual ~DarcyLMH() override;
248 
249 
250 protected:
251  /**
252  * Returns true is the fields involved in the time term have values that makes the time term zero.
253  * For time_global==true, it returns true if there are no field descriptors in the input list, so the
254  * fields )of the time ter) have their default values for whole simulation.
255  * If time_global==false (default), only the actual values are considered.
256  */
257  virtual bool zero_time_term(bool time_global=false);
258 
259  /// Solve method common to zero_time_step and update solution.
260  void solve_nonlinear();
261 
262  /**
263  * Create and preallocate MH linear system (including matrix, rhs and solution vectors)
264  */
266 
267  /**
268  * Read initial condition into solution vector.
269  * Must be called after create_linear_system.
270  *
271  * For the LMH scheme we have to be able to save edge pressures in order to
272  * restart simulation or use results of one simulation as initial condition for other one.
273  */
274 // void read_initial_condition();
275 
276  /**
277  * In some circumstances, the intial condition must be processed.
278  * It is called at the end of @p read_initial_condition().
279  * This is used in Richards equation due the update of water content.
280  */
281 // virtual void initial_condition_postprocess();
282 
283  /**
284  * Allocates linear system matrix for MH.
285  * TODO:
286  * - use general preallocation methods in DofHandler
287  */
288  void allocate_mh_matrix();
289 
290  /**
291  * Assembles linear system matrix for MH.
292  * Element by element assembly is done using dim-template assembly class.
293  * Assembles only steady part of the equation.
294  * TODO:
295  * - include time term - DONE
296  * - add support for Robin type sources
297  * - support for nonlinear solvers - assembly either residual vector, matrix, or both (using FADBAD++)
298  */
299 // void assembly_mh_matrix(MultidimAssembly& assembler);
300 
301 // void reconstruct_solution_from_schur(MultidimAssembly& assembler);
302 
303  /**
304  * Assembly or update whole linear system.
305  */
306  virtual void assembly_linear_system();
307 
308 // void set_mesh_data_for_bddc(LinSys_BDDC * bddc_ls);
309  /**
310  * Return a norm of residual vector.
311  * TODO: Introduce Equation::compute_residual() updating
312  * residual field, standard part of EqData.
313  */
314 // virtual double solution_precision() const;
315 
316  /// Print darcy flow matrix in matlab format into a file.
317  void print_matlab_matrix(string matlab_file);
318 
319  /// Get vector of all DOF indices of given component (0..side, 1..element, 2..edge)
320  std::vector<int> get_component_indices_vec(unsigned int component) const;
321 
322  /// Getter for the linear system of the 2. Schur complement.
324  { return *(eq_data_->lin_sys_schur); }
325 
326  /// Create and initialize assembly objects
327  virtual void initialize_asm();
328 
329  /// Call assemble of read_init_cond_assembly_
330  virtual void read_init_cond_asm();
331 
332  std::shared_ptr<Balance> balance_;
333 
335 
336  int size; // global size of MH matrix
337 
339 
340  // Setting of the nonlinear solver. TODO: Move to the solver class later on.
341  double tolerance_;
342  unsigned int min_n_it_;
343  unsigned int max_n_it_;
344 
345  std::shared_ptr<EqFields> eq_fields_;
346  std::shared_ptr<EqData> eq_data_;
347 
348 
349  friend class DarcyFlowMHOutput;
350  //friend class P0_CouplingAssembler;
351  //friend class P1_CouplingAssembler;
352 
353  /// general assembly objects, hold assembly objects of appropriate dimension
357 private:
358  /// Registrar of class to factory
359  static const int registrar;
360 };
361 
362 #endif //DARCY_FLOW_LMH_HH
363 //-----------------------------------------------------------------------------
364 // vim: set cindent:
365 
bc_field.hh
DarcyLMH::EqData::save_local_system_
std::vector< bool > save_local_system_
Flag for saving the local system. Currently used only in case of seepage BC.
Definition: darcy_flow_lmh.hh:203
ReadInitCondAssemblyLMH
Definition: assembly_lmh.hh:46
DarcyLMH::print_matlab_matrix
void print_matlab_matrix(string matlab_file)
Print darcy flow matrix in matlab format into a file.
Definition: darcy_flow_lmh.cc:1105
DarcyLMH::zero_time_term
virtual bool zero_time_term(bool time_global=false)
Definition: darcy_flow_lmh.cc:583
DarcyLMH::EqData::reset
void reset()
Reset data members.
Definition: darcy_flow_lmh.cc:179
LocalConstraint
Definition: local_constraint.hh:34
LinSys
Definition: la_linsys_new.hh:169
vector_mpi.hh
DarcyLMH::output_object
DarcyFlowMHOutput * output_object
Definition: darcy_flow_lmh.hh:334
time_governor.hh
Basic time management class.
LocalSystem
Definition: local_system.hh:46
DarcyLMH::EqData
Definition: darcy_flow_lmh.hh:158
GenericAssemblyBase
Definition: generic_assembly.hh:52
DarcyLMH::EqData::p_edge_solution
VectorMPI p_edge_solution
Definition: darcy_flow_lmh.hh:176
Input
Abstract linear system class.
Definition: balance.hh:40
DarcyLMH::accept_time_step
virtual void accept_time_step()
postprocess velocity field (add sources)
Definition: darcy_flow_lmh.cc:686
string.h
DarcyLMH::EqData::bc_fluxes_reconstruted
std::vector< bool > bc_fluxes_reconstruted
Flag indicating whether the fluxes for seepage BC has been reconstructed already.
Definition: darcy_flow_lmh.hh:204
Balance
Definition: balance.hh:119
DarcyLMH::create_linear_system
void create_linear_system(Input::AbstractRecord rec)
Definition: darcy_flow_lmh.cc:925
DarcyLMH::EqData::seepage_bc_systems
std::map< LongIdx, LocalSystem > seepage_bc_systems
Definition: darcy_flow_lmh.hh:180
DarcyLMH::TYPEDEF_ERR_INFO
TYPEDEF_ERR_INFO(EI_Reason, string)
DarcyLMH::init_eq_data
void init_eq_data()
Definition: darcy_flow_lmh.cc:275
DarcyLMH::EqFields::EqFields
EqFields()
Definition: darcy_flow_lmh.cc:157
field_set.hh
DarcyLMH::EqData::schur_offset_
std::array< unsigned int, 3 > schur_offset_
Index offset in the local system for the Schur complement (of dim = 1,2,3).
Definition: darcy_flow_lmh.hh:205
DarcyLMH::max_n_it_
unsigned int max_n_it_
Definition: darcy_flow_lmh.hh:343
std::vector< LocalSystem >
DarcyLMH::data_changed_
bool data_changed_
Definition: darcy_flow_lmh.hh:338
DarcyLMH::EqData::lin_sys_schur
std::shared_ptr< LinSys > lin_sys_schur
Definition: darcy_flow_lmh.hh:175
DarcyLMH::solve_time_step
void solve_time_step(bool output=true)
Solve the problem without moving to next time and without output.
Definition: darcy_flow_lmh.cc:529
DarcyLMH::update_solution
void update_solution() override
Definition: darcy_flow_lmh.cc:515
DarcyLMH::EqFields
Definition: darcy_flow_lmh.hh:151
DarcyLMH::get_input_type
static const Input::Type::Record & get_input_type()
Definition: darcy_flow_lmh.cc:104
type_base.hh
DarcyLMH::initialize_asm
virtual void initialize_asm()
Create and initialize assembly objects.
Definition: darcy_flow_lmh.cc:1379
DarcyLMH::postprocess
virtual void postprocess()
Definition: darcy_flow_lmh.cc:1037
DarcyLMH::DECLARE_EXCEPTION
DECLARE_EXCEPTION(ExcSolverDiverge,<< "Diverged nonlinear solver. Reason: "<< EI_Reason::val)
DarcyLMH::assembly_linear_system
virtual void assembly_linear_system()
Definition: darcy_flow_lmh.cc:1065
DarcyLMH::set_extra_source
void set_extra_source(const Field< 3, FieldValue< 3 >::Scalar > &extra_src)
Sets external source field (coupling with other equation).
Definition: darcy_flow_lmh.hh:244
DarcyLMH::mh_matrix_assembly_
GenericAssemblyBase * mh_matrix_assembly_
Definition: darcy_flow_lmh.hh:355
DarcyLMH::min_n_it_
unsigned int min_n_it_
Definition: darcy_flow_lmh.hh:342
exceptions.hh
Element
Definition: elements.h:40
DarcyLMH::type_field_descriptor
static const Input::Type::Record & type_field_descriptor()
Definition: darcy_flow_lmh.cc:89
DarcyMH::EqData
Definition: darcy_flow_mh.hh:204
DarcyLMH::zero_time_step
void zero_time_step() override
Definition: darcy_flow_lmh.cc:472
EquationBase::mesh
Mesh & mesh()
Definition: equation.hh:178
DarcyLMH::EqData::loc_edge_dofs
std::array< std::vector< unsigned int >, 3 > loc_edge_dofs
Definition: darcy_flow_lmh.hh:198
DarcyLMH::eq_fields
EqFields & eq_fields()
Definition: darcy_flow_lmh.hh:236
Intersection
Definition: intersectionquadrature.hh:146
Input::Record
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
DarcyLMH::eq_data
EqData & eq_data()
Definition: darcy_flow_lmh.hh:237
DarcyLMH::EqData::postprocess_solution_
std::vector< arma::vec > postprocess_solution_
Definition: darcy_flow_lmh.hh:196
type_generic.hh
DarcyLMH::EqData::dh_p_
std::shared_ptr< SubDOFHandlerMultiDim > dh_p_
DOF handler represents DOFs of element pressure.
Definition: darcy_flow_lmh.hh:166
DarcyLMH::solve_nonlinear
void solve_nonlinear()
Solve method common to zero_time_step and update solution.
Definition: darcy_flow_lmh.cc:592
DarcyLMH::EqData::loc_system_
std::vector< LocalSystem > loc_system_
Following data members are stored in vectors, one item for every cell.
Definition: darcy_flow_lmh.hh:194
DarcyMH::EqFields
Definition: darcy_flow_mh.hh:143
DarcyLMH::balance_
std::shared_ptr< Balance > balance_
Definition: darcy_flow_lmh.hh:332
DarcyLMH::registrar
static const int registrar
Registrar of class to factory.
Definition: darcy_flow_lmh.hh:359
DarcyLMH::set_extra_storativity
void set_extra_storativity(const Field< 3, FieldValue< 3 >::Scalar > &extra_stor)
Sets external storarivity field (coupling with other equation).
Definition: darcy_flow_lmh.hh:240
TimeGovernor
Basic time management functionality for unsteady (and steady) solvers (class Equation).
Definition: time_governor.hh:317
DarcyLMH::EqData::EqData
EqData()
Definition: darcy_flow_lmh.cc:163
DarcyLMH::solved_time
virtual double solved_time() override
Definition: darcy_flow_lmh.cc:252
field_values.hh
Input::AbstractRecord
Accessor to the polymorphic input data of a type given by an AbstracRecord object.
Definition: accessors.hh:458
input_exception.hh
mesh.h
Input::Type::Selection
Template for classes storing finite set of named values.
Definition: type_selection.hh:65
std::map< LongIdx, LocalSystem >
DarcyLMH::initialize_specific
virtual void initialize_specific()
Definition: darcy_flow_lmh.cc:429
DarcyLMH::eq_data_
std::shared_ptr< EqData > eq_data_
Definition: darcy_flow_lmh.hh:346
DarcyLMH::~DarcyLMH
virtual ~DarcyLMH() override
Definition: darcy_flow_lmh.cc:1344
DarcyLMH::EqData::nonlinear_iteration_
unsigned int nonlinear_iteration_
Definition: darcy_flow_lmh.hh:185
DarcyLMH::initialize
void initialize() override
Definition: darcy_flow_lmh.cc:318
Input::Type::Record
Record type proxy class.
Definition: type_record.hh:182
DarcyLMH::get_component_indices_vec
std::vector< int > get_component_indices_vec(unsigned int component) const
Get vector of all DOF indices of given component (0..side, 1..element, 2..edge)
Definition: darcy_flow_lmh.cc:1365
DarcyLMH::EqData::init
void init()
Initialize vectors, ...
Definition: darcy_flow_lmh.cc:169
DarcyLMH::EqData::time_step_
double time_step_
Definition: darcy_flow_lmh.hh:173
DarcyLMH::lin_sys_schur
LinSys & lin_sys_schur()
Getter for the linear system of the 2. Schur complement.
Definition: darcy_flow_lmh.hh:323
DarcyLMH::EqData::use_steady_assembly_
bool use_steady_assembly_
Definition: darcy_flow_lmh.hh:170
Mesh
Definition: mesh.h:362
DarcyLMH::tolerance_
double tolerance_
Definition: darcy_flow_lmh.hh:341
DarcyLMH::EqData::loc_side_dofs
std::array< std::vector< unsigned int >, 3 > loc_side_dofs
Definition: darcy_flow_lmh.hh:197
DarcyLMH::get_mh_mortar_selection
static const Input::Type::Selection & get_mh_mortar_selection()
Selection for enum MortarMethod.
Definition: darcy_flow_lmh.cc:81
DarcyLMH::allocate_mh_matrix
void allocate_mh_matrix()
Definition: darcy_flow_lmh.cc:750
DarcyLMH
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
Definition: darcy_flow_lmh.hh:133
darcy_flow_mh.hh
mixed-hybrid model of linear Darcy flow, possibly unsteady.
DarcyLMH::reconstruct_schur_assembly_
GenericAssemblyBase * reconstruct_schur_assembly_
Definition: darcy_flow_lmh.hh:356
DarcyLMH::size
int size
Definition: darcy_flow_lmh.hh:336
DarcyLMH::read_init_cond_assembly_
GenericAssembly< ReadInitCondAssemblyLMH > * read_init_cond_assembly_
general assembly objects, hold assembly objects of appropriate dimension
Definition: darcy_flow_lmh.hh:354
DarcyLMH::EqData::p_edge_solution_previous
VectorMPI p_edge_solution_previous
Definition: darcy_flow_lmh.hh:177
DarcyLMH::EqData::loc_constraint_
std::vector< LocalConstraint > loc_constraint_
Definition: darcy_flow_lmh.hh:195
DarcyLMH::EqData::balance_
std::shared_ptr< Balance > balance_
Shared Balance object.
Definition: darcy_flow_lmh.hh:183
DarcyLMH::DECLARE_INPUT_EXCEPTION
DECLARE_INPUT_EXCEPTION(ExcMissingTimeGovernor,<< "Missing the key 'time', obligatory for the transient problems.")
VectorMPI
Definition: vector_mpi.hh:43
Field
Class template representing a field with values dependent on: point, element, and region.
Definition: field.hh:91
DarcyLMH::read_init_cond_asm
virtual void read_init_cond_asm()
Call assemble of read_init_cond_assembly_.
Definition: darcy_flow_lmh.cc:1386
DarcyFlowMHOutput
Definition: darcy_flow_mh_output.hh:74
DarcyLMH::eq_fields_
std::shared_ptr< EqFields > eq_fields_
Definition: darcy_flow_lmh.hh:345
DarcyLMH::EqData::loc_ele_dof
std::array< unsigned int, 3 > loc_ele_dof
Definition: darcy_flow_lmh.hh:199
DarcyLMH::output_data
virtual void output_data() override
Write computed fields.
Definition: darcy_flow_lmh.cc:694
GenericAssembly
Generic class of assemblation.
Definition: generic_assembly.hh:160
DarcyLMH::DarcyLMH
DarcyLMH(Mesh &mesh, const Input::Record in_rec, TimeGovernor *tm=nullptr)
CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL.
Definition: darcy_flow_lmh.cc:202
darcy_flow_interface.hh
field.hh
DarcyFlowInterface
Definition: darcy_flow_interface.hh:15
FieldValue
Definition: field_values.hh:645
DarcyLMH::EqData::p_edge_solution_previous_time
VectorMPI p_edge_solution_previous_time
Definition: darcy_flow_lmh.hh:178
IntDim
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19