![]()  | 
  
    LifeV
    
   | 
 
implements a mixed-hybrid FE Darcy solver for transient problems More...
#include <DarcySolverTransient.hpp>
 Inheritance diagram for DarcySolverTransient< MeshType >:
 Collaboration diagram for DarcySolverTransient< MeshType >:Constructors & destructor | |
| DarcySolverTransient () | |
| Constructor for the class.  More... | |
| virtual | ~DarcySolverTransient () | 
| Virtual destructor.  More... | |
methods | |
| virtual void | setup () | 
| Set up the linear solver and the preconditioner for the linear system.  More... | |
| virtual void | solve () | 
| Solve the problem.  More... | |
Update and set methods | |
| void | setInitialPrimal (const scalarFctPtr_Type &primalInitialFct) | 
| Initialize primal solution.  More... | |
| void | setMass (const scalarFctPtr_Type &massFct) | 
| Set mass matrix.  More... | |
Protected methods | |
| virtual void | localMatrixComputation (const UInt &iElem, MatrixElemental &elmatMix, MatrixElemental &elmatReactionTerm) | 
| Compute element matrices.  More... | |
| virtual void | localVectorComputation (const UInt &iElem, VectorElemental &elvecMix) | 
| Computes local vectors.  More... | |
| virtual void | postComputePrimalAndDual () | 
| Do some computation after the calculation of the primal and dual variable.  More... | |
| void | setupTime () | 
| Setup the time data.  More... | |
Time advance stuff | |
| timeAdvancePtr_Type | M_timeAdvance | 
| Time advance.  More... | |
| vectorPtr_Type | M_rhsTimeAdvance | 
| Right hand side coming from the time advance scheme.  More... | |
| matrixElementalContainer_Type | M_localMassMatrix | 
| Local mass matrices.  More... | |
Private Constructors | |
| DarcySolverTransient (const darcySolverTransient_Type &) | |
| Inhibited copy constructor.  More... | |
Private Operators | |
| darcySolverTransient_Type & | operator= (const darcySolverTransient_Type &) | 
| Inhibited assign operator.  More... | |
Data of the problem | |
| scalarFctPtr_Type | M_primalFieldInitialFct | 
| Initial time primal variable.  More... | |
| scalarFctPtr_Type | M_massFct | 
| Mass function, it does not depend on time.  More... | |
Algebraic stuff | |
| bool | M_reusePrec | 
| Boolean that indicates if the preconditioner is re-used or not.  More... | |
| bool | M_updated | 
| Boolean that indicates if the matrix is updated for the current iteration.  More... | |
| UInt | M_maxIterSolver | 
| Interger storing the max number of solver iteration with preconditioner recomputing.  More... | |
| bool | M_recomputeMatrix | 
| Boolean that indicates if the matrix is recomputed for the current iteration.  More... | |
Additional Inherited Members | |
  Public Types inherited from DarcySolverLinear< MeshType > | |
| typedef MeshType | mesh_Type | 
| Typedef for mesh template.  More... | |
| typedef LinearSolver | solver_Type | 
| Typedef for solver template.  More... | |
| typedef DarcySolverLinear< mesh_Type > | darcySolver_Type | 
| Self typedef.  More... | |
| typedef DarcyData< mesh_Type > | data_Type | 
| Typedef for the data type.  More... | |
| typedef std::shared_ptr< data_Type > | dataPtr_Type | 
| Shared pointer for the data type.  More... | |
| typedef BCHandler | bcHandler_Type | 
| Boundary condition handler.  More... | |
| typedef std::shared_ptr< bcHandler_Type > | bcHandlerPtr_Type | 
| Shared pointer to a boundary condition handler.  More... | |
| typedef std::shared_ptr< Epetra_Comm > | commPtr_Type | 
| Shared pointer to a MPI communicator.  More... | |
| typedef MapEpetra | map_Type | 
| Map type.  More... | |
| typedef std::shared_ptr< Displayer > | displayerPtr_Type | 
| Shared pointer to a displayer.  More... | |
| typedef FESpace< mesh_Type, map_Type > | fESpace_Type | 
| Finite element space.  More... | |
| typedef std::shared_ptr< fESpace_Type > | fESpacePtr_Type | 
| Shared pointer to a finite element space.  More... | |
| typedef FEScalarField< mesh_Type, map_Type > | scalarField_Type | 
| Scalar field.  More... | |
| typedef std::shared_ptr< scalarField_Type > | scalarFieldPtr_Type | 
| Shared pointer to a scalar field.  More... | |
| typedef FEVectorField< mesh_Type, map_Type > | vectorField_Type | 
| Vector field.  More... | |
| typedef std::shared_ptr< vectorField_Type > | vectorFieldPtr_Type | 
| Shared pointer to a scalar field.  More... | |
| typedef FEFunction< mesh_Type, map_Type, Real > | scalarFct_Type | 
| Scalar value function.  More... | |
| typedef std::shared_ptr< scalarFct_Type > | scalarFctPtr_Type | 
| Shared pointer to a scalar value function.  More... | |
| typedef FEFunction< mesh_Type, map_Type, Vector > | vectorFct_Type | 
| Vector value function.  More... | |
| typedef std::shared_ptr< vectorFct_Type > | vectorFctPtr_Type | 
| Shared pointer to a vector value function.  More... | |
| typedef FEFunction< mesh_Type, map_Type, Matrix > | matrixFct_Type | 
| Matrix value funcion.  More... | |
| typedef std::shared_ptr< matrixFct_Type > | matrixFctPtr_Type | 
| Shared pointer to a matrix value function.  More... | |
| typedef solver_Type::matrix_Type | matrix_Type | 
| Sparse and distributed matrix.  More... | |
| typedef solver_Type::matrixPtr_Type | matrixPtr_Type | 
| Shared pointer to a sparse and distributed matrix.  More... | |
| typedef solver_Type::vector_Type | vector_Type | 
| Distributed vector.  More... | |
| typedef solver_Type::vectorPtr_Type | vectorPtr_Type | 
| Shared pointer to a distributed vector.  More... | |
| typedef solver_Type::preconditionerPtr_Type | preconditionerPtr_Type | 
| Shared pointer to the preconditioner.  More... | |
  Public Member Functions inherited from DarcySolverLinear< MeshType > | |
| DarcySolverLinear () | |
| Constructor for the class.  More... | |
| virtual | ~DarcySolverLinear () | 
| Virtual destructor.  More... | |
| void | buildSystem () | 
| Build the global hybrid system, the right hand and apply the boundary conditions.  More... | |
| void | solveLinearSystem () | 
| Solve the global hybrid system.  More... | |
| void | computePrimalAndDual () | 
| Compute primal and dual variables from the hybrid variable as a post process.  More... | |
| void | setData (dataPtr_Type &data) | 
| Set the data.  More... | |
| void | setBoundaryConditions (bcHandlerPtr_Type &bcHandler) | 
| Set the boundary conditions.  More... | |
| void | setScalarSource (const scalarFctPtr_Type &scalarSourceFct) | 
| Set scalar source term.  More... | |
| void | setVectorSource (const vectorFctPtr_Type &vectorSourceFct) | 
| Set vector source term.  More... | |
| virtual void | setInversePermeability (const matrixFctPtr_Type &invPermFct) | 
| Set inverse diffusion tensor.  More... | |
| void | setReactionTerm (const scalarFctPtr_Type &reactionTermFct) | 
| Set the coefficient for the reaction term.  More... | |
| void | setHybridField (const scalarFieldPtr_Type &hybridField) | 
| Set the hybrid field vector.  More... | |
| void | setPrimalField (const scalarFieldPtr_Type &primalField) | 
| Set the primal field vector.  More... | |
| void | setDualField (const scalarFieldPtr_Type &dualField) | 
| Set the dual field vector.  More... | |
| void | setFields (const scalarFieldPtr_Type &dualField, const scalarFieldPtr_Type &primalField, const scalarFieldPtr_Type &hybridField) | 
| Set the fields.  More... | |
| void | setDisplayer (const displayerPtr_Type &displayer) | 
| Set the displayer and, implicitly, the communicator.  More... | |
| void | setCommunicator (const commPtr_Type &comm) | 
| Set the communicator and, implicitly, the displayer.  More... | |
| const scalarFieldPtr_Type & | hybridFieldPtr () const | 
| Returns pointer to the hybrid solution field.  More... | |
| scalarFieldPtr_Type & | hybridFieldPtr () | 
| Returns pointer to the hybrid solution field.  More... | |
| const scalarFieldPtr_Type & | primalFieldPtr () const | 
| Returns pointer to the primal solution field.  More... | |
| scalarFieldPtr_Type & | primalFieldPtr () | 
| Returns pointer to the primal solution field.  More... | |
| const scalarFieldPtr_Type & | dualFieldPtr () const | 
| Returns pointer to the dual solution field.  More... | |
| scalarFieldPtr_Type & | dualFieldPtr () | 
| Returns pointer to the dual solution field.  More... | |
| const vectorPtr_Type & | residualPtr () const | 
| vectorPtr_Type & | residualPtr () | 
| const bcHandlerPtr_Type & | boundaryConditionHandlerPtr () const | 
| Returns boundary conditions handler.  More... | |
| bcHandlerPtr_Type & | boundaryConditionHandlerPtr () | 
| Returns boundary conditions handler.  More... | |
| const map_Type & | getMap () const | 
| Returns Epetra local map.  More... | |
| const commPtr_Type & | getCommPtr () const | 
| Returns Epetra communicator.  More... | |
  Protected Member Functions inherited from DarcySolverLinear< MeshType > | |
| DarcySolverLinear (const darcySolver_Type &) | |
| Inhibited copy constructor.  More... | |
| darcySolver_Type & | operator= (const darcySolver_Type &) | 
| Inhibited assign operator.  More... | |
| virtual void | preLoopElementsComputation () | 
| Perform all the operations before doing the loop on volume elements.  More... | |
| void | computeConstantMatrices (MatrixElemental &elmatMix) | 
| Pre-computes local (element independant) matrices.  More... | |
| void | staticCondensation (MatrixElemental &localMatrixHybrid, VectorElemental &localVectorHybrid, MatrixElemental &elmatMix, MatrixElemental &elmatReactionTerm, VectorElemental &elvecMix) | 
| Performs static condensation.  More... | |
| void | localComputePrimalAndDual (VectorElemental &localSolution, MatrixElemental &elmatMix, MatrixElemental &elmatReactionTerm, VectorElemental &elvecMix) | 
| Compute locally, as a post process, the primal and dual variable given the hybrid.  More... | |
| virtual void | resetVariables () | 
| Update all problem variables.  More... | |
| void | applyBoundaryConditions () | 
| Apply the boundary condition.  More... | |
| template<typename MatrixType > | |
| void | symmetrizeMatrix (Int N, MatrixType &A) | 
| Make matrix symmetric.  More... | |
  Protected Attributes inherited from DarcySolverLinear< MeshType > | |
| displayerPtr_Type | M_displayer | 
| Displayer for parallel cout.  More... | |
| dataPtr_Type | M_data | 
| Data for Darcy solvers.  More... | |
| scalarFctPtr_Type | M_scalarSourceFct | 
| Source function.  More... | |
| vectorFctPtr_Type | M_vectorSourceFct | 
| Vector source function.  More... | |
| scalarFctPtr_Type | M_reactionTermFct | 
| Reaction term function.  More... | |
| matrixFctPtr_Type | M_inversePermeabilityFct | 
| Inverse of the permeability tensor.  More... | |
| bcHandlerPtr_Type | M_boundaryConditionHandler | 
| Bondary conditions handler.  More... | |
| scalarFieldPtr_Type | M_primalField | 
| Primal solution.  More... | |
| scalarFieldPtr_Type | M_dualField | 
| Dual solution.  More... | |
| scalarFieldPtr_Type | M_hybridField | 
| Hybrid solution.  More... | |
| matrixPtr_Type | M_matrHybrid | 
| Hybrid matrix.  More... | |
| vectorPtr_Type | M_rhs | 
| Right hand side.  More... | |
| vectorPtr_Type | M_residual | 
| Residual.  More... | |
| solver_Type | M_linearSolver | 
| Linear solver.  More... | |
| preconditionerPtr_Type | M_prec | 
| Epetra preconditioner for the linear system.  More... | |
implements a mixed-hybrid FE Darcy solver for transient problems
This class implements a transient Darcy solver in the finite time interval 
. 
 The classical time dependent Darcy formulation is a couple of differential equations of first order with the unknowns 
, being the pressure or the primal unknown, and 
, being the Darcy velocity or the flux or the dual unknown, such that 
The data in the system are:
 the primal variable at initial time;  
 the mass term, which does not depend on time;  
 the permeability tensor;  
 the scalar source term;  
 the vector source term;  
 the reaction coefficient;  
 the subset of the boundary of 
 with Dirichlet boundary conditions;  
 Dirichlet boundary condition data;  
 that is the part of the boundary of 
 with Robin, or Neumann, boundary conditions;  
 and 
 Robin boundary condition data;  We suppose that 
 and 
. 
 Using the hybridization procedure, and introducing a new variable, we may split the problem from the entire domain to problems in each elements of the triangulation 
, then we may write the weak formulation for the dual problem. The new variable is the hybrid unknown 
, e.g. "the trace" of pressure of the primal unknown, it is the Lagrange multipliers forcing the continuity of the flux across each faces or edges in the triangulation. We introduce the following functional spaces, the first is the space of the primal variable, the third for the dual variable, the fourth for the hybrid variable. 
 Introducing the following bilinear forms and functionals at fixed time 
 
 we obtain the Darcy problem in the weak form: find 
 such that 
 At semi-discrete level. i.e. only space discretization is performed, we introduce the polynomial space, of degree 
, that approximate the finite dimensional spaces introduced above 
, 
 and 
 
 where 
 is the space of polynomial of degree 
 in the element 
, 
 is the space of polynomial of Raviart-Thomas of degrees 
 in the element 
 and 
 is the space of polynomial of degree 
 definite on each of the boundary of the element 
 and discontinuous from one edge to the other. 
 The finite dimensional problem is: find 
 such that 
 To obatin a fully discrete problem we discretize the time derivative via BDF (backward differentiation formulae) schemes of order 
, with 
 as time step. Choosing 
 we obtain the following system for each 
, with 
 
 To solve the problem we use the static condensation procedure, i.e. the unknowns in the discrete weak system are not independent and 
, 
 may be written in function of 
 alone. We introduce the following local matrices 
 where we avoid to write the dependence on the triangle 
 and on the current time step 
 in all the matrices and vectors. 
 bre local matrix formulation of the finite dimensional problem is 
Or alternatively
Introducing the local hybrid matrix and local hybrid right hand side
Imposing that at each edge or face the hybrid unknown is single value we obtain a linear system for the hybrid unknown
We recover the primal and dual variable as a post-process from the hybrid variable at element level, so we have
 and the vector 
, because all the boundary conditions are imposed via BCHandler class. Definition at line 251 of file DarcySolverTransient.hpp.
Typedef for mesh template.
Definition at line 261 of file DarcySolverTransient.hpp.
Self typedef.
Definition at line 264 of file DarcySolverTransient.hpp.
| typedef DarcySolverLinear< mesh_Type > darcySolverLinear_Type | 
Darcy solver class.
Definition at line 267 of file DarcySolverTransient.hpp.
Typedef for the data type.
Definition at line 270 of file DarcySolverTransient.hpp.
Shared pointer for the data type.
Definition at line 273 of file DarcySolverTransient.hpp.
Shared pointer to a scalar value function.
Definition at line 276 of file DarcySolverTransient.hpp.
Shared pointer to a distributed vector.
Definition at line 279 of file DarcySolverTransient.hpp.
Distributed vector.
Definition at line 282 of file DarcySolverTransient.hpp.
Shared pointer to a scalar field.
Definition at line 285 of file DarcySolverTransient.hpp.
Scalar field.
Definition at line 288 of file DarcySolverTransient.hpp.
| typedef std::vector< MatrixElemental > matrixElementalContainer_Type | 
Container of matrix elemental.
Definition at line 291 of file DarcySolverTransient.hpp.
| typedef TimeAdvanceBDF< vector_Type > timeAdvance_Type | 
Time advance scheme.
Definition at line 294 of file DarcySolverTransient.hpp.
| typedef std::shared_ptr< timeAdvance_Type > timeAdvancePtr_Type | 
Shared pointer to a time advance scheme.
Definition at line 297 of file DarcySolverTransient.hpp.
Constructor for the class.
Definition at line 467 of file DarcySolverTransient.hpp.
 Here is the caller graph for this function:
      
  | 
  inlinevirtual | 
Virtual destructor.
Definition at line 308 of file DarcySolverTransient.hpp.
      
  | 
  private | 
Inhibited copy constructor.
      
  | 
  virtual | 
Set up the linear solver and the preconditioner for the linear system.
Reimplemented from DarcySolverLinear< MeshType >.
Reimplemented in DarcySolverTransientNonLinear< MeshType >.
Definition at line 488 of file DarcySolverTransient.hpp.
      
  | 
  virtual | 
Solve the problem.
Reimplemented from DarcySolverLinear< MeshType >.
Reimplemented in DarcySolverTransientNonLinear< MeshType >.
Definition at line 552 of file DarcySolverTransient.hpp.
| void setInitialPrimal | ( | const scalarFctPtr_Type & | primalInitialFct | ) | 
Initialize primal solution.
Set the initial value function for the primal variable and compute the primal variable at step zero.
| primalInitial | The primal initial function. | 
Definition at line 526 of file DarcySolverTransient.hpp.
| void setMass | ( | const scalarFctPtr_Type & | massFct | ) | 
Set mass matrix.
Set the mass term, the default source term is the function one. By defaul it does not depend on time.
| mass | Mass term for the problem. | 
Loop on all the volume elements.
Definition at line 573 of file DarcySolverTransient.hpp.
      
  | 
  protectedvirtual | 
Compute element matrices.
Call the Darcy solver localMatrixComputation method and compute the mass matrix for the time dependent term.
| iElem | Id of the current geometrical element. | 
| elmatMix | The local matrix in mixed form. | 
| elmatReactionTerm | The local matrix for the reaction term. | 
Reimplemented from DarcySolverLinear< MeshType >.
Reimplemented in DarcySolverTransientNonLinear< MeshType >.
Definition at line 624 of file DarcySolverTransient.hpp.
      
  | 
  protectedvirtual | 
Computes local vectors.
Call the Darc_y solver localVectorComputation method and compute the additional scalar vector for the time dependent term.
| iElem | Id of the current geometrical element. | 
| elvecMix | The local vector in mixed form. | 
Reimplemented from DarcySolverLinear< MeshType >.
Reimplemented in DarcySolverTransientNonLinear< MeshType >.
Definition at line 650 of file DarcySolverTransient.hpp.
      
  | 
  inlineprotectedvirtual | 
Do some computation after the calculation of the primal and dual variable.
Save into the time advance scheme the computed primal solution.
Reimplemented from DarcySolverLinear< MeshType >.
Definition at line 378 of file DarcySolverTransient.hpp.
      
  | 
  protected | 
Setup the time data.
Definition at line 503 of file DarcySolverTransient.hpp.
 Here is the caller graph for this function:
      
  | 
  private | 
Inhibited assign operator.
      
  | 
  protected | 
Time advance.
Definition at line 394 of file DarcySolverTransient.hpp.
      
  | 
  protected | 
Right hand side coming from the time advance scheme.
Definition at line 397 of file DarcySolverTransient.hpp.
      
  | 
  private | 
Initial time primal variable.
Definition at line 426 of file DarcySolverTransient.hpp.
      
  | 
  private | 
Mass function, it does not depend on time.
Definition at line 429 of file DarcySolverTransient.hpp.
      
  | 
  private | 
Local mass matrices.
Definition at line 438 of file DarcySolverTransient.hpp.
      
  | 
  private | 
Boolean that indicates if the preconditioner is re-used or not.
Definition at line 447 of file DarcySolverTransient.hpp.
      
  | 
  private | 
Boolean that indicates if the matrix is updated for the current iteration.
Definition at line 450 of file DarcySolverTransient.hpp.
      
  | 
  private | 
Interger storing the max number of solver iteration with preconditioner recomputing.
Definition at line 453 of file DarcySolverTransient.hpp.
      
  | 
  private | 
Boolean that indicates if the matrix is recomputed for the current iteration.
Definition at line 456 of file DarcySolverTransient.hpp.