LifeV
|
implements a mixed-hybrid FE Darcy solver. More...
#include <DarcySolverLinear.hpp>
Public Types | |
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... | |
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_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... | |
Private Constructors | |
DarcySolverLinear (const darcySolver_Type &) | |
Inhibited copy constructor. More... | |
Private Operators | |
darcySolver_Type & | operator= (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... | |
implements a mixed-hybrid FE Darcy solver.
This class implements a Darcy solver.
The classical formulation of this problem is a couple of differential equations of first order with two unknowns: , being the pressure or the primal unknown, and , being the Darcy velocity or the total flux or the dual unknown, such that
The first equation is the Darcy equation and the second equation is the conservation law.
The data in the system are:
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
we obtain the Darcy problem in the weak form: find such that
At discrete level 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 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 in all the matrices and vectors.
The 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
Definition at line 245 of file DarcySolverLinear.hpp.
Typedef for mesh template.
Definition at line 254 of file DarcySolverLinear.hpp.
typedef LinearSolver solver_Type |
Typedef for solver template.
Definition at line 257 of file DarcySolverLinear.hpp.
typedef DarcySolverLinear< mesh_Type > darcySolver_Type |
Self typedef.
Definition at line 260 of file DarcySolverLinear.hpp.
Typedef for the data type.
Definition at line 263 of file DarcySolverLinear.hpp.
typedef std::shared_ptr< data_Type > dataPtr_Type |
Shared pointer for the data type.
Definition at line 266 of file DarcySolverLinear.hpp.
typedef BCHandler bcHandler_Type |
Boundary condition handler.
Definition at line 269 of file DarcySolverLinear.hpp.
typedef std::shared_ptr< bcHandler_Type > bcHandlerPtr_Type |
Shared pointer to a boundary condition handler.
Definition at line 272 of file DarcySolverLinear.hpp.
typedef std::shared_ptr< Epetra_Comm > commPtr_Type |
Shared pointer to a MPI communicator.
Definition at line 275 of file DarcySolverLinear.hpp.
Map type.
Definition at line 278 of file DarcySolverLinear.hpp.
typedef std::shared_ptr< Displayer > displayerPtr_Type |
Shared pointer to a displayer.
Definition at line 281 of file DarcySolverLinear.hpp.
typedef FESpace< mesh_Type, map_Type > fESpace_Type |
Finite element space.
Definition at line 284 of file DarcySolverLinear.hpp.
typedef std::shared_ptr< fESpace_Type > fESpacePtr_Type |
Shared pointer to a finite element space.
Definition at line 287 of file DarcySolverLinear.hpp.
typedef FEScalarField< mesh_Type, map_Type > scalarField_Type |
Scalar field.
Definition at line 290 of file DarcySolverLinear.hpp.
typedef std::shared_ptr< scalarField_Type > scalarFieldPtr_Type |
Shared pointer to a scalar field.
Definition at line 293 of file DarcySolverLinear.hpp.
typedef FEVectorField< mesh_Type, map_Type > vectorField_Type |
Vector field.
Definition at line 296 of file DarcySolverLinear.hpp.
typedef std::shared_ptr< vectorField_Type > vectorFieldPtr_Type |
Shared pointer to a scalar field.
Definition at line 299 of file DarcySolverLinear.hpp.
typedef FEFunction< mesh_Type, map_Type, Real > scalarFct_Type |
Scalar value function.
Definition at line 302 of file DarcySolverLinear.hpp.
typedef std::shared_ptr< scalarFct_Type > scalarFctPtr_Type |
Shared pointer to a scalar value function.
Definition at line 305 of file DarcySolverLinear.hpp.
typedef FEFunction< mesh_Type, map_Type, Vector > vectorFct_Type |
Vector value function.
Definition at line 308 of file DarcySolverLinear.hpp.
typedef std::shared_ptr< vectorFct_Type > vectorFctPtr_Type |
Shared pointer to a vector value function.
Definition at line 311 of file DarcySolverLinear.hpp.
typedef FEFunction< mesh_Type, map_Type, Matrix > matrixFct_Type |
Matrix value funcion.
Definition at line 314 of file DarcySolverLinear.hpp.
typedef std::shared_ptr< matrixFct_Type > matrixFctPtr_Type |
Shared pointer to a matrix value function.
Definition at line 317 of file DarcySolverLinear.hpp.
typedef solver_Type::matrix_Type matrix_Type |
Sparse and distributed matrix.
Definition at line 320 of file DarcySolverLinear.hpp.
Shared pointer to a sparse and distributed matrix.
Definition at line 323 of file DarcySolverLinear.hpp.
typedef solver_Type::vector_Type vector_Type |
Distributed vector.
Definition at line 326 of file DarcySolverLinear.hpp.
Shared pointer to a distributed vector.
Definition at line 329 of file DarcySolverLinear.hpp.
Shared pointer to the preconditioner.
Definition at line 332 of file DarcySolverLinear.hpp.
|
inline |
Constructor for the class.
Definition at line 340 of file DarcySolverLinear.hpp.
|
inlinevirtual |
Virtual destructor.
Definition at line 343 of file DarcySolverLinear.hpp.
|
protected |
Inhibited copy constructor.
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.
|
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.
void solveLinearSystem | ( | ) |
Solve the global hybrid system.
Definition at line 954 of file DarcySolverLinear.hpp.
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.
|
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.
|
inline |
Set the data.
data | Data for the problem. |
Definition at line 374 of file DarcySolverLinear.hpp.
|
inline |
Set the boundary conditions.
bcHandler | Boundary condition handler for the problem. |
Definition at line 383 of file DarcySolverLinear.hpp.
|
inline |
Set scalar source term.
scalarSourceFct | Vector source term for the problem. |
Definition at line 392 of file DarcySolverLinear.hpp.
|
inline |
Set vector source term.
vectorSourceFct | Vector source term for the problem. |
Definition at line 401 of file DarcySolverLinear.hpp.
|
inlinevirtual |
Set inverse diffusion tensor.
Set the inverse of diffusion tensor, the default inverse of permeability is the identity matrix.
invPermFct | Inverse of the permeability tensor for the problem. |
Reimplemented in DarcySolverNonLinear< MeshType >, and DarcySolverTransientNonLinear< MeshType >.
Definition at line 411 of file DarcySolverLinear.hpp.
|
inline |
Set the coefficient for the reaction term.
reactionTermFct | Reaction term for the problem. |
Definition at line 421 of file DarcySolverLinear.hpp.
|
inline |
Set the hybrid field vector.
hybrid | Constant scalarFieldPtr_Type reference of the hybrid vector. |
Definition at line 430 of file DarcySolverLinear.hpp.
|
inline |
Set the primal field vector.
primalField | Constant scalarFieldPtr_Type reference of the primal solution. |
Definition at line 439 of file DarcySolverLinear.hpp.
|
inline |
Set the dual field vector.
dualField | Constant scalarFieldPtr_Type reference of the dual solution. |
Definition at line 448 of file DarcySolverLinear.hpp.
|
inline |
Set the fields.
dualField | Constant scalarFieldPtr_Type reference of the dual solution. |
primalField | Constant scalarFieldPtr_Type reference of the primal solution. |
hybrid | Constant scalarFieldPtr_Type reference of the hybrid vector. |
Definition at line 459 of file DarcySolverLinear.hpp.
|
inline |
Set the displayer and, implicitly, the communicator.
displayer | Constant displayerPtr_Type reference of the displayer. |
Definition at line 477 of file DarcySolverLinear.hpp.
|
inline |
Set the communicator and, implicitly, the displayer.
comm | Constant commPtr_Type reference of the communicator. |
Definition at line 486 of file DarcySolverLinear.hpp.
|
inline |
Returns pointer to the hybrid solution field.
Definition at line 500 of file DarcySolverLinear.hpp.
|
inline |
Returns pointer to the hybrid solution field.
Definition at line 509 of file DarcySolverLinear.hpp.
|
inline |
Returns pointer to the primal solution field.
Definition at line 518 of file DarcySolverLinear.hpp.
|
inline |
Returns pointer to the primal solution field.
Definition at line 527 of file DarcySolverLinear.hpp.
|
inline |
Returns pointer to the dual solution field.
Definition at line 536 of file DarcySolverLinear.hpp.
|
inline |
Returns pointer to the dual solution field.
Definition at line 545 of file DarcySolverLinear.hpp.
|
inline |
Returns the pointer of the residual vector.
Definition at line 554 of file DarcySolverLinear.hpp.
|
inline |
Returns the pointer of the residual vector.
Definition at line 563 of file DarcySolverLinear.hpp.
|
inline |
Returns boundary conditions handler.
Definition at line 572 of file DarcySolverLinear.hpp.
|
inline |
Returns boundary conditions handler.
Definition at line 581 of file DarcySolverLinear.hpp.
|
inline |
Returns Epetra local map.
Definition at line 590 of file DarcySolverLinear.hpp.
|
inline |
Returns Epetra communicator.
Definition at line 599 of file DarcySolverLinear.hpp.
|
protected |
Inhibited assign operator.
|
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.
|
protected |
Pre-computes local (element independant) matrices.
Compute all the local matrices that are independant from the geometrical element.
elmatMix | The local matrix in mixed form. |
Definition at line 1125 of file DarcySolverLinear.hpp.
|
protectedvirtual |
Computes local (element dependent) matrices.
Locally update the current finite element for the dual finite element space, then compute the Hdiv mass matrix.
iElem | Id of the current geometrical element. |
elmatMix | The local matrix in mixed form. |
elmatReactionTerm | The local matrix for the reaction term. |
Reimplemented in DarcySolverTransientNonLinear< MeshType >, and DarcySolverTransient< MeshType >.
Definition at line 1174 of file DarcySolverLinear.hpp.
|
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.
iElem | Id of the current geometrical element. |
elvecMix | The local vector in mixed form. |
Reimplemented in DarcySolverTransientNonLinear< MeshType >, and DarcySolverTransient< MeshType >.
Definition at line 1225 of file DarcySolverLinear.hpp.
|
protected |
Performs static condensation.
Locally eliminate pressure and velocity DOFs, create the local hybrid matrix and local hybrid right hand side.
localMatrixHybrid | The matrix which will store the hybrid local matrix. |
localVectorHybrid | The vector which will store the hybrid local vector. |
elmatMix | The local matrix in mixed form. |
elmatReactionTerm | The local matrix for the reaction term. |
elvecMix | The 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.
|
protected |
Compute locally, as a post process, the primal and dual variable given the hybrid.
localSolution | A vector which stores the dual, primal and hybrid local solution. |
elmatMix | The local matrix in mixed form. |
elmatReactionTerm | The local matrix for the reaction term. |
elvecMix | The 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.
|
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.
|
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.
|
protected |
Apply the boundary condition.
Definition at line 1659 of file DarcySolverLinear.hpp.
|
protected |
Make matrix symmetric.
Transform a symmetric matrix that is stored only in the lower triangular part in a full symmetric matrix.
N | The size of the matrix A. |
A | The matrix to be reordered. |
Definition at line 1693 of file DarcySolverLinear.hpp.
|
protected |
Displayer for parallel cout.
Definition at line 734 of file DarcySolverLinear.hpp.
|
protected |
Data for Darcy solvers.
Definition at line 743 of file DarcySolverLinear.hpp.
|
protected |
Source function.
Definition at line 746 of file DarcySolverLinear.hpp.
|
protected |
Vector source function.
Definition at line 749 of file DarcySolverLinear.hpp.
|
protected |
Reaction term function.
Definition at line 752 of file DarcySolverLinear.hpp.
|
protected |
Inverse of the permeability tensor.
Definition at line 755 of file DarcySolverLinear.hpp.
|
protected |
Bondary conditions handler.
Definition at line 758 of file DarcySolverLinear.hpp.
|
protected |
Primal solution.
Definition at line 767 of file DarcySolverLinear.hpp.
|
protected |
Dual solution.
Definition at line 770 of file DarcySolverLinear.hpp.
|
protected |
Hybrid solution.
Definition at line 773 of file DarcySolverLinear.hpp.
|
protected |
Hybrid matrix.
Definition at line 782 of file DarcySolverLinear.hpp.
|
protected |
Right hand side.
Definition at line 785 of file DarcySolverLinear.hpp.
|
protected |
Residual.
Definition at line 788 of file DarcySolverLinear.hpp.
|
protected |
Linear solver.
Definition at line 791 of file DarcySolverLinear.hpp.
|
protected |
Epetra preconditioner for the linear system.
Definition at line 794 of file DarcySolverLinear.hpp.