MFEM v4.7.0
Finite element discretization library
|
Bramble-Pasciak Solver for Darcy equation. More...
#include <bramble_pasciak.hpp>
Public Member Functions | |
BramblePasciakSolver (ParBilinearForm *mVarf, ParMixedBilinearForm *bVarf, const BPSParameters ¶m) | |
System and mass preconditioner are constructed from bilinear forms. | |
BramblePasciakSolver (HypreParMatrix &M, HypreParMatrix &B, HypreParMatrix &Q, Solver &M0, Solver &M1, const BPSParameters ¶m) | |
System and mass preconditioner are user-provided. | |
virtual void | Mult (const Vector &x, Vector &y) const |
Operator application: y=A(x) . | |
virtual void | SetOperator (const Operator &op) |
Set/update the solver for the given operator. | |
void | SetEssZeroDofs (const Array< int > &dofs) |
virtual int | GetNumIterations () const |
Public Member Functions inherited from mfem::blocksolvers::DarcySolver | |
DarcySolver (int size0, int size1) | |
Public Member Functions inherited from mfem::Solver | |
Solver (int s=0, bool iter_mode=false) | |
Initialize a square Solver with size s. | |
Solver (int h, int w, bool iter_mode=false) | |
Initialize a Solver with height h and width w. | |
Public Member Functions inherited from mfem::Operator | |
void | InitTVectors (const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const |
Initializes memory for true vectors of linear system. | |
Operator (int s=0) | |
Construct a square Operator with given size s (default 0). | |
Operator (int h, int w) | |
Construct an Operator with the given height (output size) and width (input size). | |
int | Height () const |
Get the height (size of output) of the Operator. Synonym with NumRows(). | |
int | NumRows () const |
Get the number of rows (size of output) of the Operator. Synonym with Height(). | |
int | Width () const |
Get the width (size of input) of the Operator. Synonym with NumCols(). | |
int | NumCols () const |
Get the number of columns (size of input) of the Operator. Synonym with Width(). | |
virtual MemoryClass | GetMemoryClass () const |
Return the MemoryClass preferred by the Operator. | |
virtual void | MultTranspose (const Vector &x, Vector &y) const |
Action of the transpose operator: y=A^t(x) . The default behavior in class Operator is to generate an error. | |
virtual void | AddMult (const Vector &x, Vector &y, const real_t a=1.0) const |
Operator application: y+=A(x) (default) or y+=a*A(x) . | |
virtual void | AddMultTranspose (const Vector &x, Vector &y, const real_t a=1.0) const |
Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x) . | |
virtual void | ArrayMult (const Array< const Vector * > &X, Array< Vector * > &Y) const |
Operator application on a matrix: Y=A(X) . | |
virtual void | ArrayMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y) const |
Action of the transpose operator on a matrix: Y=A^t(X) . | |
virtual void | ArrayAddMult (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const |
Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X) . | |
virtual void | ArrayAddMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const |
Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X) . | |
virtual Operator & | GetGradient (const Vector &x) const |
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error. | |
virtual void | AssembleDiagonal (Vector &diag) const |
Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operators. In some cases, only an approximation of the diagonal is computed. | |
virtual const Operator * | GetProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. | |
virtual const Operator * | GetRestriction () const |
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. | |
virtual const Operator * | GetOutputProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. | |
virtual const Operator * | GetOutputRestrictionTranspose () const |
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators. | |
virtual const Operator * | GetOutputRestriction () const |
Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity. | |
void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0) |
Form a constrained linear system using a matrix-free approach. | |
void | FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B) |
Form a column-constrained linear system using a matrix-free approach. | |
virtual void | RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x) |
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear system obtained from Operator::FormLinearSystem() or Operator::FormRectangularLinearSystem(). | |
void | FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A) |
Return in A a parallel (on truedofs) version of this square operator. | |
void | FormRectangularSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A) |
Return in A a parallel (on truedofs) version of this rectangular operator (including constraints). | |
void | FormDiscreteOperator (Operator *&A) |
Return in A a parallel (on truedofs) version of this rectangular operator. | |
void | PrintMatlab (std::ostream &out, int n, int m=0) const |
Prints operator with input size n and output size m in Matlab format. | |
virtual void | PrintMatlab (std::ostream &out) const |
Prints operator in Matlab format. | |
virtual | ~Operator () |
Virtual destructor. | |
Type | GetType () const |
Return the type ID of the Operator class. | |
Static Public Member Functions | |
static HypreParMatrix * | ConstructMassPreconditioner (ParBilinearForm &mVarf, real_t alpha=0.5) |
Assemble a preconditioner for the mass matrix. | |
Additional Inherited Members | |
Public Types inherited from mfem::Operator | |
enum | DiagonalPolicy { DIAG_ZERO , DIAG_ONE , DIAG_KEEP } |
Defines operator diagonal policy upon elimination of rows and/or columns. More... | |
enum | Type { ANY_TYPE , MFEM_SPARSEMAT , Hypre_ParCSR , PETSC_MATAIJ , PETSC_MATIS , PETSC_MATSHELL , PETSC_MATNEST , PETSC_MATHYPRE , PETSC_MATGENERIC , Complex_Operator , MFEM_ComplexSparseMat , Complex_Hypre_ParCSR , Complex_DenseMat , MFEM_Block_Matrix , MFEM_Block_Operator } |
Enumeration defining IDs for some classes derived from Operator. More... | |
Public Attributes inherited from mfem::Solver | |
bool | iterative_mode |
If true, use the second argument of Mult() as an initial guess. | |
Protected Member Functions inherited from mfem::Operator | |
void | FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout) |
see FormSystemOperator() | |
void | FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout) |
see FormRectangularSystemOperator() | |
Operator * | SetupRAP (const Operator *Pi, const Operator *Po) |
Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt". | |
Protected Attributes inherited from mfem::blocksolvers::DarcySolver | |
Array< int > | offsets_ |
Protected Attributes inherited from mfem::Operator | |
int | height |
Dimension of the output / number of rows in the matrix. | |
int | width |
Dimension of the input / number of columns in the matrix. | |
Bramble-Pasciak Solver for Darcy equation.
Bramble-Pasciak Solver for Darcy equation.
The basic idea is to precondition the mass matrix M with a s.p.d. matrix Q such that M - Q remains s.p.d. Then we can transform the block operator into a s.p.d. operator under a modified inner product. In particular, this enable us to implement modified versions of CG iterations, that rely on efficient applications of the required transformations.
We offer a mass preconditioner based on a rescaling of the diagonal of the element mass matrices M_T.
We consider Q_T := alpha * lambda_min * D_T, where D_T := diag(M_T), and lambda_min is the smallest eigenvalue of the following problem
M_T x = lambda * D_T x.
Alpha is a parameter that is strictly between 0 and 1.
For more details, see:
Definition at line 122 of file bramble_pasciak.hpp.
mfem::blocksolvers::BramblePasciakSolver::BramblePasciakSolver | ( | ParBilinearForm * | mVarf, |
ParMixedBilinearForm * | bVarf, | ||
const BPSParameters & | param ) |
System and mass preconditioner are constructed from bilinear forms.
Bramble-Pasciak Solver.
Definition at line 22 of file bramble_pasciak.cpp.
mfem::blocksolvers::BramblePasciakSolver::BramblePasciakSolver | ( | HypreParMatrix & | M, |
HypreParMatrix & | B, | ||
HypreParMatrix & | Q, | ||
Solver & | M0, | ||
Solver & | M1, | ||
const BPSParameters & | param ) |
System and mass preconditioner are user-provided.
Definition at line 46 of file bramble_pasciak.cpp.
|
static |
Assemble a preconditioner for the mass matrix.
Mass preconditioner corresponds to a local re-scaling based on the smallest eigenvalue of the generalized eigenvalue problem locally on each element T: M_T x_T = lambda_T diag(M_T) x_T. We set Q_T = alpha * min(lambda_T) * diag(M_T), 0 < alpha < 1.
Definition at line 135 of file bramble_pasciak.cpp.
|
inlinevirtual |
Implements mfem::blocksolvers::DarcySolver.
Definition at line 167 of file bramble_pasciak.hpp.
Operator application: y=A(x)
.
Implements mfem::Operator.
Definition at line 207 of file bramble_pasciak.cpp.
|
inline |
Definition at line 166 of file bramble_pasciak.hpp.
|
inlinevirtual |
Set/update the solver for the given operator.
Implements mfem::Solver.
Definition at line 165 of file bramble_pasciak.hpp.