LifeV
DarcySolverLinear< MeshType > Class Template Reference

implements a mixed-hybrid FE Darcy solver. More...

#include <DarcySolverLinear.hpp>

+ Inheritance diagram for DarcySolverLinear< MeshType >:
+ Collaboration diagram for DarcySolverLinear< MeshType >:

Public Types

typedef MeshType mesh_Type
 Typedef for mesh template. More...
 
typedef LinearSolver solver_Type
 Typedef for solver template. More...
 
typedef DarcySolverLinear< mesh_TypedarcySolver_Type
 Self typedef. More...
 
typedef DarcyData< mesh_Typedata_Type
 Typedef for the data type. More...
 
typedef std::shared_ptr< data_TypedataPtr_Type
 Shared pointer for the data type. More...
 
typedef BCHandler bcHandler_Type
 Boundary condition handler. More...
 
typedef std::shared_ptr< bcHandler_TypebcHandlerPtr_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< DisplayerdisplayerPtr_Type
 Shared pointer to a displayer. More...
 
typedef FESpace< mesh_Type, map_TypefESpace_Type
 Finite element space. More...
 
typedef std::shared_ptr< fESpace_TypefESpacePtr_Type
 Shared pointer to a finite element space. More...
 
typedef FEScalarField< mesh_Type, map_TypescalarField_Type
 Scalar field. More...
 
typedef std::shared_ptr< scalarField_TypescalarFieldPtr_Type
 Shared pointer to a scalar field. More...
 
typedef FEVectorField< mesh_Type, map_TypevectorField_Type
 Vector field. More...
 
typedef std::shared_ptr< vectorField_TypevectorFieldPtr_Type
 Shared pointer to a scalar field. More...
 
typedef FEFunction< mesh_Type, map_Type, RealscalarFct_Type
 Scalar value function. More...
 
typedef std::shared_ptr< scalarFct_TypescalarFctPtr_Type
 Shared pointer to a scalar value function. More...
 
typedef FEFunction< mesh_Type, map_Type, VectorvectorFct_Type
 Vector value function. More...
 
typedef std::shared_ptr< vectorFct_TypevectorFctPtr_Type
 Shared pointer to a vector value function. More...
 
typedef FEFunction< mesh_Type, map_Type, MatrixmatrixFct_Type
 Matrix value funcion. More...
 
typedef std::shared_ptr< matrixFct_TypematrixFctPtr_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...
 

Constructors and destructor

 DarcySolverLinear ()
 Constructor for the class. More...
 
virtual ~DarcySolverLinear ()
 Virtual destructor. More...
 

Methods

void buildSystem ()
 Build the global hybrid system, the right hand and apply the boundary conditions. More...
 
virtual void setup ()
 Set up the linear solver and the preconditioner for the linear system. 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...
 
virtual void solve ()
 Solve the Darcy problem grouping other public methods. More...
 

Set Methods

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...
 

Get Methods

const scalarFieldPtr_TypehybridFieldPtr () const
 Returns pointer to the hybrid solution field. More...
 
scalarFieldPtr_TypehybridFieldPtr ()
 Returns pointer to the hybrid solution field. More...
 
const scalarFieldPtr_TypeprimalFieldPtr () const
 Returns pointer to the primal solution field. More...
 
scalarFieldPtr_TypeprimalFieldPtr ()
 Returns pointer to the primal solution field. More...
 
const scalarFieldPtr_TypedualFieldPtr () const
 Returns pointer to the dual solution field. More...
 
scalarFieldPtr_TypedualFieldPtr ()
 Returns pointer to the dual solution field. More...
 
const vectorPtr_TyperesidualPtr () const
 
vectorPtr_TyperesidualPtr ()
 
const bcHandlerPtr_TypeboundaryConditionHandlerPtr () const
 Returns boundary conditions handler. More...
 
bcHandlerPtr_TypeboundaryConditionHandlerPtr ()
 Returns boundary conditions handler. More...
 
const map_TypegetMap () const
 Returns Epetra local map. More...
 
const commPtr_TypegetCommPtr () const
 Returns Epetra communicator. More...
 

Private Constructors

 DarcySolverLinear (const darcySolver_Type &)
 Inhibited copy constructor. More...
 

Private Operators

darcySolver_Typeoperator= (const darcySolver_Type &)
 Inhibited assign operator. More...
 

Protected Methods

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...
 
virtual void localMatrixComputation (const UInt &iElem, MatrixElemental &elmatMix, MatrixElemental &elmatReactionTerm)
 Computes local (element dependent) matrices. More...
 
virtual void localVectorComputation (const UInt &iElem, VectorElemental &elvecMix)
 Computes local (element dependent) vectors. 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 postComputePrimalAndDual ()
 Do some computation after the calculation of the primal and dual variable. 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...
 

Parallel stuff

displayerPtr_Type M_displayer
 Displayer for parallel cout. More...
 

Data of the problem

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...
 

Solution fields stuff

scalarFieldPtr_Type M_primalField
 Primal solution. More...
 
scalarFieldPtr_Type M_dualField
 Dual solution. More...
 
scalarFieldPtr_Type M_hybridField
 Hybrid solution. More...
 

Algebraic stuff

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...
 

Detailed Description

template<typename MeshType>
class LifeV::DarcySolverLinear< MeshType >

implements a mixed-hybrid FE Darcy solver.

Author
A. Fumagalli aless.nosp@m.io.f.nosp@m.umaga.nosp@m.lli@.nosp@m.mail..nosp@m.poli.nosp@m.mi.it
See also
For applications related to two-phase flow see [8]

This class implements a Darcy solver.

The classical formulation of this problem is a couple of differential equations of first order with two unknowns: $ p \in C^1 (\Omega ) $, being the pressure or the primal unknown, and $ \sigma \in (C^1( \Omega ) )^n $, being the Darcy velocity or the total flux or the dual unknown, such that

\[ \left\{ \begin{array}{l l l} \Lambda^{-1} \sigma + \nabla p = f_v & \mathrm{in} & \Omega\,, \vspace{0.2cm} \\ \nabla \cdot \sigma + \pi p - f = 0 & \mathrm{in} & \Omega\,, \vspace{0.2cm} \\ p = g_D & \mathrm{on} & \Gamma_D\,,\vspace{0.2cm} \\ \sigma \cdot n + h p = g_R & \mathrm{on} & \Gamma_R\,. \end{array} \right. \]

The first equation is the Darcy equation and the second equation is the conservation law.
The data in the system are:

  • $ \Lambda $ the permeability tensor;
  • $ f $ the scalar source term;
  • $ f_v $ the vector source term;
  • $ \pi $ the reaction coefficient;
  • $ \Gamma_D $ the subset of the boundary of $ \Omega $ with Dirichlet boundary conditions;
  • $ g_D $ Dirichlet boundary condition data;
  • $ \Gamma_R $ that is the part of the boundary of $ \Omega $ with Robin, or Neumann, boundary conditions;
  • $ h $ and $ g_R $ Robin boundary condition data;

We suppose that $ \partial \Omega = \Gamma_D \cup \Gamma_R $ and $ \Gamma_D \cap \Gamma_R = \emptyset $.
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 $ \mathcal{T}_h $, then we may write the weak formulation for the dual problem. The new variable is the hybrid unknown $ \lambda $, 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.

\[ \begin{array}{l} V = L^2 (\Omega ) \,,\vspace{0.2cm} \\ H(div, K) = \displaystyle \left\{ \tau \in ( L^2(K))^n : \nabla \cdot \tau \in L^2(K)\right\}\,, \vspace{0.2cm}\\ Z = \displaystyle \left\{ \tau \in L^2(\Omega) : \tau\vert_K \in H(div, K) \, \forall K \in \mathcal{T}_h \right\}\,, \vspace{0.2cm}\\ \Lambda = \displaystyle \left\{ \lambda \in \prod_{K \in \mathcal{T}_h} H^{1/2} (\partial K): \lambda_K = \lambda_{K'} \ ,\, \mathrm{on} \,\, e_{K\cap K'} \, \forall K \in \mathcal{T}_h,\, \lambda = g_D \,\, \mathrm{on} \,\, \Gamma_D \right\}\,. \end{array} \]

Introducing the following bilinear forms and functionals

\[ \begin{array}{l l} a(\sigma, \tau) = \displaystyle \sum_{K \in \mathcal{T}_h}\int_K \Lambda^{-1} \sigma \cdot \tau \,, & b(p, \tau) = \displaystyle -\sum_{K \in \mathcal{T}_h} \int_K p \nabla \cdot \tau\,, \vspace{0.2cm}\\ c(\lambda, \tau) = \displaystyle \sum_{K \in \mathcal{T}_h} \int_{\partial K} \lambda \tau \cdot n\,,& d(p, q) = \displaystyle \sum_{K \in \mathcal{T}_h} \int_{\partial K} \pi p q \,,\vspace{0.2cm}\\ h(\lambda, \mu) = \displaystyle \sum_{K \in \mathcal{T}_h} \int_{\partial K \cap \Gamma_R} h \mu \lambda \,,& F(q) = \displaystyle \sum_{K \in \mathcal{T}_h} \int_K f q\,,\vspace{0.2cm}\\ F_v(\tau) = \displaystyle \sum_{K \in \mathcal{T}_h} \int_K f_v \tau \,, \vspace{0.2cm}\,,& G(\mu) = \displaystyle \sum_{K \in \mathcal{T}_h} \int_{\partial K \cap \Gamma_R} g \mu\,, \end{array} \]

we obtain the Darcy problem in the weak form: find $ (\sigma, \, p, \, \lambda) \in Z \times V \times \Lambda $ such that

\[ \left\{ \begin{array}{l l} a(\sigma, \tau) + b(p, \tau) + c(\lambda, \tau) = F_v(\tau)\,, & \forall \tau \in Z \,,\vspace{0.2cm}\\ b(q, \sigma) - d(p,q) = - F(q)\,, & \forall q \in V \,,\vspace{0.2cm}\\ c(\mu, \sigma) + h(\lambda, \mu) = G(\mu) \,, & \forall \mu \in \Lambda\,. \end{array} \right. \]

At discrete level we introduce the polynomial space, of degree $ r $, that approximate the finite dimensional spaces introduced above $ V_h \subset V $, $ Z_h \subset Z $ and $\Lambda_h \subset \Lambda $

\[ \begin{array}{l} V_h = \displaystyle \left\{ v_h \in V: v_h|_K \in P_r (K)\, \forall K \in \mathcal{T}_h \right\}\,, \vspace{0.2cm}\\ Z_h = \displaystyle \left\{ \tau_h \in Z: \tau_h|_K \in RT_r(K) \, \forall K \in \mathcal{T}_h \right\} \,, \vspace{0.2cm} \\ \Lambda_h = \displaystyle \left\{ \lambda_h \in \Lambda: \lambda_h|_{\partial K} \in R_r( \partial K ) \, \forall K \in \mathcal{T}_h\right\}\,, \end{array} \]

where $ P_r(K) $ is the space of polynomial of degree $ r $ in the element $ K $, $ RT_r(K) $ is the space of polynomial of Raviart-Thomas of degrees $ r $ in the element $ K $ and $ R_r(\partial K) $ is the space of polynomial of degree $ r $ definite on each of the boundary of the element $ K $ and discontinuous from one edge to the other.
The finite dimensional problem is: find $ (\sigma_h,\, p_h, \, \lambda_h) \in Z_h \times V_h \times \Lambda_h $ such that

\[ \left\{ \begin{array}{l l} a(\sigma_h, \tau_h) + b(p_h, \tau_h) + c(\lambda_h, \tau_h) = F_v(\tau_h) \,, & \forall \tau_h \in Z_h \,,\vspace{0.2cm}\\ b(q_h, \sigma_h) - d(p_h, q_h)= -F(q_h)\,, & \forall q_h \in V_h \,,\vspace{0.2cm}\\ c(\mu_h, \sigma_h) + h(\lambda_h, \mu_h) = G(\mu_h) \,, & \forall \mu_h \in \Lambda_h\,. \end{array} \right. \]

To solve the problem we use the static condensation procedure, i.e. the unknowns in the discrete weak system are not independent and $ p_K $, $\sigma_K $ may be written in function of $ \lambda_K $ alone. We introduce the following local matrices

\[ \begin{array}{l l l} \left[ A \right]_{ij} = \displaystyle \int_K \Lambda^{-1} \psi_j \cdot \psi_i \,, & \left[ B \right]_{ij} = \displaystyle - \int_K \phi_j \nabla \cdot \psi_i \,, & \left[ C \right]_{ij} = \displaystyle \int_{\partial K} \xi_i \psi_j \cdot n \,, \vspace{0.2cm} \\ \left[ D \right]_{ij} = \displaystyle \int_K \pi \phi_j \phi_i\,, & \left[ H \right]_{ij} = \displaystyle \int_{\partial K \cap \Gamma_R} h \xi_i \xi_j\,, & \left[ F \right]_{j} = \displaystyle \int_K f \phi_j\,, \vspace{0.2cm} \\ \left[ F_v \right]_{j} = \displaystyle \int_K f_v \cdot \psi_j\,,& \left[ G \right]_{j} = \displaystyle \int_{\partial K \cap \Gamma_R } g \xi_j\,, \end{array} \]

where we avoid to write the dependence on the triangle $ K $ in all the matrices and vectors.
The local matrix formulation of the finite dimensional problem is

\[ \left\{ \begin{array}{l} A \sigma_K + B p_K + C \lambda_K = F_v\,, \vspace{0.2cm} \\ B^T \sigma_K - D p_K = -F \,, \vspace{0.2cm}\\ C^T \sigma_K + H \lambda_K = G\,. \end{array} \right. \]

Or alternatively

\[ \begin{array}{l l l} \left[ \begin{array}{c c c} A & B & C \vspace{0.2cm} \\ B^T & -D & 0 \vspace{0.2cm} \\ C^T & 0 & H \end{array} \right] \, \cdot & \left[ \begin{array}{c} \sigma_K \vspace{0.2cm}\\ p_K \vspace{0.2cm}\\ \lambda_K \end{array} \right] \, & = \left[ \begin{array}{c} F_v \vspace{0.2cm}\\ -F \vspace{0.2cm}\\ G \end{array} \right]\,. \end{array} \]

Introducing the local hybrid matrix and local hybrid right hand side

\[ \begin{array}{l} L_K = - C^T A^{-1} C + C^T A^{-1} B \left( B^T A^{-1} B + D \right)^{-1} B^T A^{-1} C + H \,, \vspace{0.2cm} \\ r_K = G + C^T A^{-1} B \left( B^T A^{-1} B + D \right)^{-1} F + C^T A^{-1} B \left( B^T A^{-1} B + D \right)^{-1} B^T A^{-1} F_v - C^T A^{-1} F_v \,, \end{array} \]

Imposing that at each edge or face the hybrid unknown is single value we obtain a linear system for the hybrid unknown

\[ L \lambda = r \,. \]

We recover the primal and dual variable as a post-process from the hybrid variable at element level, so we have

\[ \begin{array}{l} p_K = \left( B^T A^{-1} B + D \right)^{-1} \left( F + B^T A^{-1} F_v - B^T A^{-1} C \lambda_K \right)\,, \vspace{0.2cm} \\ \sigma_K = -A^{-1} \left( B p_K - F_v + C \lambda_K \right) \,. \end{array} \]

Note
In the code we do not use the matrix $ H $ and the vector $ G $, because all the boundary conditions are imposed via BCHandler class.
Example of usage can be found in darcy_nonlinear and darcy_linear. Coupled with an hyperbolic solver in impes.
Todo:
Insert any scientific publications that use this solver.
Bug:
If the save flag for the exporter is setted to 0 the program fails.

Definition at line 245 of file DarcySolverLinear.hpp.

Member Typedef Documentation

◆ mesh_Type

Typedef for mesh template.

Definition at line 254 of file DarcySolverLinear.hpp.

◆ solver_Type

Typedef for solver template.

Definition at line 257 of file DarcySolverLinear.hpp.

◆ darcySolver_Type

Self typedef.

Definition at line 260 of file DarcySolverLinear.hpp.

◆ data_Type

Typedef for the data type.

Definition at line 263 of file DarcySolverLinear.hpp.

◆ dataPtr_Type

typedef std::shared_ptr< data_Type > dataPtr_Type

Shared pointer for the data type.

Definition at line 266 of file DarcySolverLinear.hpp.

◆ bcHandler_Type

Boundary condition handler.

Definition at line 269 of file DarcySolverLinear.hpp.

◆ bcHandlerPtr_Type

typedef std::shared_ptr< bcHandler_Type > bcHandlerPtr_Type

Shared pointer to a boundary condition handler.

Definition at line 272 of file DarcySolverLinear.hpp.

◆ commPtr_Type

typedef std::shared_ptr< Epetra_Comm > commPtr_Type

Shared pointer to a MPI communicator.

Definition at line 275 of file DarcySolverLinear.hpp.

◆ map_Type

Map type.

Definition at line 278 of file DarcySolverLinear.hpp.

◆ displayerPtr_Type

typedef std::shared_ptr< Displayer > displayerPtr_Type

Shared pointer to a displayer.

Definition at line 281 of file DarcySolverLinear.hpp.

◆ fESpace_Type

Finite element space.

Definition at line 284 of file DarcySolverLinear.hpp.

◆ fESpacePtr_Type

typedef std::shared_ptr< fESpace_Type > fESpacePtr_Type

Shared pointer to a finite element space.

Definition at line 287 of file DarcySolverLinear.hpp.

◆ scalarField_Type

Scalar field.

Definition at line 290 of file DarcySolverLinear.hpp.

◆ scalarFieldPtr_Type

typedef std::shared_ptr< scalarField_Type > scalarFieldPtr_Type

Shared pointer to a scalar field.

Definition at line 293 of file DarcySolverLinear.hpp.

◆ vectorField_Type

Vector field.

Definition at line 296 of file DarcySolverLinear.hpp.

◆ vectorFieldPtr_Type

typedef std::shared_ptr< vectorField_Type > vectorFieldPtr_Type

Shared pointer to a scalar field.

Definition at line 299 of file DarcySolverLinear.hpp.

◆ scalarFct_Type

Scalar value function.

Definition at line 302 of file DarcySolverLinear.hpp.

◆ scalarFctPtr_Type

typedef std::shared_ptr< scalarFct_Type > scalarFctPtr_Type

Shared pointer to a scalar value function.

Definition at line 305 of file DarcySolverLinear.hpp.

◆ vectorFct_Type

Vector value function.

Definition at line 308 of file DarcySolverLinear.hpp.

◆ vectorFctPtr_Type

typedef std::shared_ptr< vectorFct_Type > vectorFctPtr_Type

Shared pointer to a vector value function.

Definition at line 311 of file DarcySolverLinear.hpp.

◆ matrixFct_Type

Matrix value funcion.

Definition at line 314 of file DarcySolverLinear.hpp.

◆ matrixFctPtr_Type

typedef std::shared_ptr< matrixFct_Type > matrixFctPtr_Type

Shared pointer to a matrix value function.

Definition at line 317 of file DarcySolverLinear.hpp.

◆ matrix_Type

Sparse and distributed matrix.

Definition at line 320 of file DarcySolverLinear.hpp.

◆ matrixPtr_Type

Shared pointer to a sparse and distributed matrix.

Definition at line 323 of file DarcySolverLinear.hpp.

◆ vector_Type

Distributed vector.

Definition at line 326 of file DarcySolverLinear.hpp.

◆ vectorPtr_Type

Shared pointer to a distributed vector.

Definition at line 329 of file DarcySolverLinear.hpp.

◆ preconditionerPtr_Type

Shared pointer to the preconditioner.

Definition at line 332 of file DarcySolverLinear.hpp.

Constructor & Destructor Documentation

◆ DarcySolverLinear() [1/2]

DarcySolverLinear ( )
inline

Constructor for the class.

Definition at line 340 of file DarcySolverLinear.hpp.

◆ ~DarcySolverLinear()

virtual ~DarcySolverLinear ( )
inlinevirtual

Virtual destructor.

Definition at line 343 of file DarcySolverLinear.hpp.

◆ DarcySolverLinear() [2/2]

DarcySolverLinear ( const darcySolver_Type )
protected

Inhibited copy constructor.

Member Function Documentation

◆ buildSystem()

void buildSystem ( )

Build the global hybrid system, the right hand and apply the boundary conditions.

Loop on all the volume elements.

End of loop volume operation.

Definition at line 812 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ setup()

void setup ( )
virtual

Set up the linear solver and the preconditioner for the linear system.

Reimplemented in DarcySolverTransient< MeshType >, DarcySolverTransientNonLinear< MeshType >, and DarcySolverNonLinear< MeshType >.

Definition at line 933 of file DarcySolverLinear.hpp.

◆ solveLinearSystem()

void solveLinearSystem ( )

Solve the global hybrid system.

Definition at line 954 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ computePrimalAndDual()

void computePrimalAndDual ( )

Compute primal and dual variables from the hybrid variable as a post process.

Loop on all the volume elements.

End of loop on the volume elements.

Definition at line 982 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ solve()

void solve ( )
virtual

Solve the Darcy problem grouping other public methods.

Reimplemented in DarcySolverTransient< MeshType >, DarcySolverTransientNonLinear< MeshType >, and DarcySolverNonLinear< MeshType >.

Definition at line 1104 of file DarcySolverLinear.hpp.

◆ setData()

void setData ( dataPtr_Type data)
inline

Set the data.

Parameters
dataData for the problem.

Definition at line 374 of file DarcySolverLinear.hpp.

◆ setBoundaryConditions()

void setBoundaryConditions ( bcHandlerPtr_Type bcHandler)
inline

Set the boundary conditions.

Parameters
bcHandlerBoundary condition handler for the problem.

Definition at line 383 of file DarcySolverLinear.hpp.

◆ setScalarSource()

void setScalarSource ( const scalarFctPtr_Type scalarSourceFct)
inline

Set scalar source term.

Parameters
scalarSourceFctVector source term for the problem.

Definition at line 392 of file DarcySolverLinear.hpp.

◆ setVectorSource()

void setVectorSource ( const vectorFctPtr_Type vectorSourceFct)
inline

Set vector source term.

Parameters
vectorSourceFctVector source term for the problem.

Definition at line 401 of file DarcySolverLinear.hpp.

◆ setInversePermeability()

virtual void setInversePermeability ( const matrixFctPtr_Type invPermFct)
inlinevirtual

Set inverse diffusion tensor.

Set the inverse of diffusion tensor, the default inverse of permeability is the identity matrix.

Parameters
invPermFctInverse of the permeability tensor for the problem.

Reimplemented in DarcySolverNonLinear< MeshType >, and DarcySolverTransientNonLinear< MeshType >.

Definition at line 411 of file DarcySolverLinear.hpp.

◆ setReactionTerm()

void setReactionTerm ( const scalarFctPtr_Type reactionTermFct)
inline

Set the coefficient for the reaction term.

Parameters
reactionTermFctReaction term for the problem.
Note
Useful also for time discretization term.

Definition at line 421 of file DarcySolverLinear.hpp.

◆ setHybridField()

void setHybridField ( const scalarFieldPtr_Type hybridField)
inline

Set the hybrid field vector.

Parameters
hybridConstant scalarFieldPtr_Type reference of the hybrid vector.

Definition at line 430 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ setPrimalField()

void setPrimalField ( const scalarFieldPtr_Type primalField)
inline

Set the primal field vector.

Parameters
primalFieldConstant scalarFieldPtr_Type reference of the primal solution.

Definition at line 439 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ setDualField()

void setDualField ( const scalarFieldPtr_Type dualField)
inline

Set the dual field vector.

Parameters
dualFieldConstant scalarFieldPtr_Type reference of the dual solution.

Definition at line 448 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ setFields()

void setFields ( const scalarFieldPtr_Type dualField,
const scalarFieldPtr_Type primalField,
const scalarFieldPtr_Type hybridField 
)
inline

Set the fields.

Parameters
dualFieldConstant scalarFieldPtr_Type reference of the dual solution.
primalFieldConstant scalarFieldPtr_Type reference of the primal solution.
hybridConstant scalarFieldPtr_Type reference of the hybrid vector.

Definition at line 459 of file DarcySolverLinear.hpp.

◆ setDisplayer()

void setDisplayer ( const displayerPtr_Type displayer)
inline

Set the displayer and, implicitly, the communicator.

Parameters
displayerConstant displayerPtr_Type reference of the displayer.

Definition at line 477 of file DarcySolverLinear.hpp.

◆ setCommunicator()

void setCommunicator ( const commPtr_Type comm)
inline

Set the communicator and, implicitly, the displayer.

Parameters
commConstant commPtr_Type reference of the communicator.

Definition at line 486 of file DarcySolverLinear.hpp.

◆ hybridFieldPtr() [1/2]

const scalarFieldPtr_Type& hybridFieldPtr ( ) const
inline

Returns pointer to the hybrid solution field.

Returns
Constant scalarFieldPtr_Type reference of the hybrid field.

Definition at line 500 of file DarcySolverLinear.hpp.

◆ hybridFieldPtr() [2/2]

scalarFieldPtr_Type& hybridFieldPtr ( )
inline

Returns pointer to the hybrid solution field.

Returns
scalarFieldPtr_Type reference of the hybrid solution field.

Definition at line 509 of file DarcySolverLinear.hpp.

◆ primalFieldPtr() [1/2]

const scalarFieldPtr_Type& primalFieldPtr ( ) const
inline

Returns pointer to the primal solution field.

Returns
Constant scalarFieldPtr_Type reference of the primal solution field.

Definition at line 518 of file DarcySolverLinear.hpp.

◆ primalFieldPtr() [2/2]

scalarFieldPtr_Type& primalFieldPtr ( )
inline

Returns pointer to the primal solution field.

Returns
scalarFieldPtr_Type reference of the primal solution field.

Definition at line 527 of file DarcySolverLinear.hpp.

◆ dualFieldPtr() [1/2]

const scalarFieldPtr_Type& dualFieldPtr ( ) const
inline

Returns pointer to the dual solution field.

Returns
Constant scalarFieldPtr_Type reference of the dual solution field.

Definition at line 536 of file DarcySolverLinear.hpp.

◆ dualFieldPtr() [2/2]

scalarFieldPtr_Type& dualFieldPtr ( )
inline

Returns pointer to the dual solution field.

Returns
scalarFieldPtr_Type reference of the dual solution field.

Definition at line 545 of file DarcySolverLinear.hpp.

◆ residualPtr() [1/2]

const vectorPtr_Type& residualPtr ( ) const
inline

Returns the pointer of the residual vector.

Returns
Constant vectorPtr_Type reference of the residual.

Definition at line 554 of file DarcySolverLinear.hpp.

◆ residualPtr() [2/2]

vectorPtr_Type& residualPtr ( )
inline

Returns the pointer of the residual vector.

Returns
vectorPtr_Type reference of the residual.

Definition at line 563 of file DarcySolverLinear.hpp.

◆ boundaryConditionHandlerPtr() [1/2]

const bcHandlerPtr_Type& boundaryConditionHandlerPtr ( ) const
inline

Returns boundary conditions handler.

Returns
Constant reference of boundary conditions handler.

Definition at line 572 of file DarcySolverLinear.hpp.

◆ boundaryConditionHandlerPtr() [2/2]

bcHandlerPtr_Type& boundaryConditionHandlerPtr ( )
inline

Returns boundary conditions handler.

Returns
Reference of boundary conditions handler.

Definition at line 581 of file DarcySolverLinear.hpp.

◆ getMap()

const map_Type& getMap ( ) const
inline

Returns Epetra local map.

Returns
Constant MapEpetra reference of the problem.

Definition at line 590 of file DarcySolverLinear.hpp.

◆ getCommPtr()

const commPtr_Type& getCommPtr ( ) const
inline

Returns Epetra communicator.

Returns
Reference of the shared pointer of the communicator of the problem.

Definition at line 599 of file DarcySolverLinear.hpp.

◆ operator=()

darcySolver_Type& operator= ( const darcySolver_Type )
protected

Inhibited assign operator.

◆ preLoopElementsComputation()

virtual void preLoopElementsComputation ( )
inlineprotectedvirtual

Perform all the operations before doing the loop on volume elements.

Computes the element independent matrices, clean all the fields and all the algebraich stuff. Useful for extentions in child classes.

Definition at line 632 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ computeConstantMatrices()

void computeConstantMatrices ( MatrixElemental elmatMix)
protected

Pre-computes local (element independant) matrices.

Compute all the local matrices that are independant from the geometrical element.

Parameters
elmatMixThe local matrix in mixed form.

Definition at line 1125 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ localMatrixComputation()

void localMatrixComputation ( const UInt iElem,
MatrixElemental elmatMix,
MatrixElemental elmatReactionTerm 
)
protectedvirtual

Computes local (element dependent) matrices.

Locally update the current finite element for the dual finite element space, then compute the Hdiv mass matrix.

Parameters
iElemId of the current geometrical element.
elmatMixThe local matrix in mixed form.
elmatReactionTermThe local matrix for the reaction term.

Reimplemented in DarcySolverTransientNonLinear< MeshType >, and DarcySolverTransient< MeshType >.

Definition at line 1174 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ localVectorComputation()

void localVectorComputation ( const UInt iElem,
VectorElemental elvecMix 
)
protectedvirtual

Computes local (element dependent) vectors.

Locally update the current finite element for the primal and dual finite element space, then compute the Hdiv and L2 source vectors.

Parameters
iElemId of the current geometrical element.
elvecMixThe local vector in mixed form.

Reimplemented in DarcySolverTransientNonLinear< MeshType >, and DarcySolverTransient< MeshType >.

Definition at line 1225 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ staticCondensation()

void staticCondensation ( MatrixElemental localMatrixHybrid,
VectorElemental localVectorHybrid,
MatrixElemental elmatMix,
MatrixElemental elmatReactionTerm,
VectorElemental elvecMix 
)
protected

Performs static condensation.

Locally eliminate pressure and velocity DOFs, create the local hybrid matrix and local hybrid right hand side.

Parameters
localMatrixHybridThe matrix which will store the hybrid local matrix.
localVectorHybridThe vector which will store the hybrid local vector.
elmatMixThe local matrix in mixed form.
elmatReactionTermThe local matrix for the reaction term.
elvecMixThe local vector in mixed form.

Matrix operations.

End of matrix operations.

Vector operations.

End of vector operations.

Definition at line 1268 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ localComputePrimalAndDual()

void localComputePrimalAndDual ( VectorElemental localSolution,
MatrixElemental elmatMix,
MatrixElemental elmatReactionTerm,
VectorElemental elvecMix 
)
protected

Compute locally, as a post process, the primal and dual variable given the hybrid.

Parameters
localSolutionA vector which stores the dual, primal and hybrid local solution.
elmatMixThe local matrix in mixed form.
elmatReactionTermThe local matrix for the reaction term.
elvecMixThe local vector in mixed form.

End of matrix operations.

Vector operations, computation of primal and dual variable.

Computation of the primal variable

Computation of the dual variable.

End of vector computation.

Definition at line 1451 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ postComputePrimalAndDual()

virtual void postComputePrimalAndDual ( )
inlineprotectedvirtual

Do some computation after the calculation of the primal and dual variable.

The function is empty while it is useful in derived classes, i.e. DarcySolverTransient.

Reimplemented in DarcySolverTransient< MeshType >.

Definition at line 701 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ resetVariables()

void resetVariables ( )
protectedvirtual

Update all problem variables.

Update all the variables of the problem before the construction of the global hybrid matrix, e.g. reset the global hybrid matrix. It is principally used for a time dependent derived class, in fact we want to clear each time step all the variables.

Reimplemented in DarcySolverNonLinear< MeshType >, and DarcySolverTransientNonLinear< MeshType >.

Definition at line 1630 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ applyBoundaryConditions()

void applyBoundaryConditions ( )
protected

Apply the boundary condition.

Definition at line 1659 of file DarcySolverLinear.hpp.

+ Here is the caller graph for this function:

◆ symmetrizeMatrix()

void symmetrizeMatrix ( Int  N,
MatrixType &  A 
)
protected

Make matrix symmetric.

Transform a symmetric matrix that is stored only in the lower triangular part in a full symmetric matrix.

Parameters
NThe size of the matrix A.
Template Parameters
AThe matrix to be reordered.

Definition at line 1693 of file DarcySolverLinear.hpp.

Field Documentation

◆ M_displayer

displayerPtr_Type M_displayer
protected

Displayer for parallel cout.

Definition at line 734 of file DarcySolverLinear.hpp.

◆ M_data

dataPtr_Type M_data
protected

Data for Darcy solvers.

Definition at line 743 of file DarcySolverLinear.hpp.

◆ M_scalarSourceFct

scalarFctPtr_Type M_scalarSourceFct
protected

Source function.

Definition at line 746 of file DarcySolverLinear.hpp.

◆ M_vectorSourceFct

vectorFctPtr_Type M_vectorSourceFct
protected

Vector source function.

Definition at line 749 of file DarcySolverLinear.hpp.

◆ M_reactionTermFct

scalarFctPtr_Type M_reactionTermFct
protected

Reaction term function.

Definition at line 752 of file DarcySolverLinear.hpp.

◆ M_inversePermeabilityFct

matrixFctPtr_Type M_inversePermeabilityFct
protected

Inverse of the permeability tensor.

Definition at line 755 of file DarcySolverLinear.hpp.

◆ M_boundaryConditionHandler

bcHandlerPtr_Type M_boundaryConditionHandler
protected

Bondary conditions handler.

Definition at line 758 of file DarcySolverLinear.hpp.

◆ M_primalField

scalarFieldPtr_Type M_primalField
protected

Primal solution.

Definition at line 767 of file DarcySolverLinear.hpp.

◆ M_dualField

scalarFieldPtr_Type M_dualField
protected

Dual solution.

Definition at line 770 of file DarcySolverLinear.hpp.

◆ M_hybridField

scalarFieldPtr_Type M_hybridField
protected

Hybrid solution.

Definition at line 773 of file DarcySolverLinear.hpp.

◆ M_matrHybrid

matrixPtr_Type M_matrHybrid
protected

Hybrid matrix.

Definition at line 782 of file DarcySolverLinear.hpp.

◆ M_rhs

vectorPtr_Type M_rhs
protected

Right hand side.

Definition at line 785 of file DarcySolverLinear.hpp.

◆ M_residual

vectorPtr_Type M_residual
protected

Residual.

Definition at line 788 of file DarcySolverLinear.hpp.

◆ M_linearSolver

solver_Type M_linearSolver
protected

Linear solver.

Definition at line 791 of file DarcySolverLinear.hpp.

◆ M_prec

preconditionerPtr_Type M_prec
protected

Epetra preconditioner for the linear system.

Definition at line 794 of file DarcySolverLinear.hpp.


The documentation for this class was generated from the following file: