MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
mfem::PetscParMatrix Class Reference

Wrapper for PETSc's matrix class. More...

#include <petsc.hpp>

Inheritance diagram for mfem::PetscParMatrix:
[legend]
Collaboration diagram for mfem::PetscParMatrix:
[legend]

Public Member Functions

 PetscParMatrix ()
 Create an empty matrix to be used as a reference to an existing matrix. More...
 
 PetscParMatrix (Mat a, bool ref=false)
 Creates PetscParMatrix out of PETSc's Mat. More...
 
 PetscParMatrix (const HypreParMatrix *ha, Operator::Type tid)
 Convert a HypreParMatrix ha to a PetscParMatrix in the given PETSc format tid. More...
 
 PetscParMatrix (MPI_Comm comm, const Operator *op, Operator::Type tid)
 Convert an mfem::Operator into a PetscParMatrix in the given PETSc format tid. More...
 
 PetscParMatrix (MPI_Comm comm, PetscInt glob_size, PetscInt *row_starts, SparseMatrix *diag, Operator::Type tid)
 Creates block-diagonal square parallel matrix. More...
 
 PetscParMatrix (MPI_Comm comm, PetscInt global_num_rows, PetscInt global_num_cols, PetscInt *row_starts, PetscInt *col_starts, SparseMatrix *diag, Operator::Type tid)
 Creates block-diagonal rectangular parallel matrix. More...
 
virtual ~PetscParMatrix ()
 Calls PETSc's destroy function. More...
 
void Mult (double a, const Vector &x, double b, Vector &y) const
 Matvec: y = a A x + b y. More...
 
void MultTranspose (double a, const Vector &x, double b, Vector &y) const
 Matvec transpose: y = a A^T x + b y. More...
 
virtual void Mult (const Vector &x, Vector &y) const
 Operator application: y=A(x). More...
 
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. More...
 
MPI_Comm GetComm () const
 Get the associated MPI communicator. More...
 
 operator Mat () const
 Typecasting to PETSc's Mat type. More...
 
PetscInt GetNumRows () const
 Returns the local number of rows. More...
 
PetscInt GetNumCols () const
 Returns the local number of columns. More...
 
PetscInt M () const
 Returns the global number of rows. More...
 
PetscInt N () const
 Returns the global number of columns. More...
 
PetscInt GetGlobalNumRows () const
 Returns the global number of rows. More...
 
PetscInt GetGlobalNumCols () const
 Returns the global number of columns. More...
 
PetscInt NNZ () const
 Returns the number of nonzeros. More...
 
PetscParVectorGetX () const
 Returns the inner vector in the domain of A (it creates it if needed) More...
 
PetscParVectorGetY () const
 Returns the inner vector in the range of A (it creates it if needed) More...
 
PetscParMatrixTranspose (bool action=false)
 Returns the transpose of the PetscParMatrix. More...
 
void Print (const char *fname=NULL, bool binary=false) const
 Prints the matrix (to stdout if fname is NULL) More...
 
void operator*= (double s)
 Scale all entries by s: A_scaled = s*A. More...
 
void EliminateRowsCols (const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
 Eliminate rows and columns from the matrix, and rows from the vector B. Modify B with the BC values in X. More...
 
PetscParMatrixEliminateRowsCols (const Array< int > &rows_cols)
 Eliminate rows and columns from the matrix and store the eliminated elements in a new matrix Ae (returned). More...
 
void MakeRef (const PetscParMatrix &master)
 Makes this object a reference to another PetscParMatrix. More...
 
Mat ReleaseMat (bool dereference)
 Release the PETSc Mat object. If dereference is true, decrement the refcount of the Mat object. More...
 
Type GetType () const
 
Assignment operators
PetscParMatrixoperator= (const PetscParMatrix &B)
 
PetscParMatrixoperator= (const HypreParMatrix &B)
 
PetscParMatrixoperator+= (const PetscParMatrix &B)
 
- Public Member Functions inherited from mfem::Operator
 Operator (int s=0)
 Construct a square Operator with given size s (default 0). More...
 
 Operator (int h, int w)
 Construct an Operator with the given height (output size) and width (input size). More...
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows(). More...
 
int NumRows () const
 Get the number of rows (size of output) of the Operator. Synonym with Height(). More...
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols(). More...
 
int NumCols () const
 Get the number of columns (size of input) of the Operator. Synonym with Width(). More...
 
virtual OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error. More...
 
virtual const OperatorGetProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More...
 
virtual const OperatorGetRestriction () const
 Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More...
 
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. More...
 
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(). More...
 
void PrintMatlab (std::ostream &out, int n=0, int m=0) const
 Prints operator with input size n and output size m in Matlab format. More...
 
virtual ~Operator ()
 Virtual destructor. More...
 

Protected Member Functions

void Init ()
 Initialize with defaults. Does not initialize inherited members. More...
 
void Destroy ()
 Delete all owned data. Does not perform re-initialization with defaults. More...
 
void MakeWrapper (MPI_Comm comm, const Operator *op, Mat *B)
 Creates a wrapper around a mfem::Operator op using PETSc's MATSHELL object and returns the Mat in B. More...
 
void ConvertOperator (MPI_Comm comm, const Operator &op, Mat *B, bool assembled)
 Convert an mfem::Operator into a Mat B; op can be destroyed. More...
 

Protected Attributes

Mat A
 The actual PETSc object. More...
 
PetscParVectorX
 Auxiliary vectors for typecasting. More...
 
PetscParVectorY
 
- Protected Attributes inherited from mfem::Operator
int height
 Dimension of the output / number of rows in the matrix. More...
 
int width
 Dimension of the input / number of columns in the matrix. More...
 

Friends

class PetscLinearSolver
 
class PetscPreconditioner
 

Additional Inherited Members

- Public Types inherited from mfem::Operator
enum  Type {
  MFEM_SPARSEMAT, HYPRE_PARCSR, PETSC_MATAIJ, PETSC_MATIS,
  PETSC_MATSHELL, PETSC_MATNEST
}
 Enumeration defining IDs for some classes derived from Operator. More...
 

Detailed Description

Wrapper for PETSc's matrix class.

Definition at line 128 of file petsc.hpp.

Constructor & Destructor Documentation

mfem::PetscParMatrix::PetscParMatrix ( )

Create an empty matrix to be used as a reference to an existing matrix.

Definition at line 366 of file petsc.cpp.

mfem::PetscParMatrix::PetscParMatrix ( Mat  a,
bool  ref = false 
)

Creates PetscParMatrix out of PETSc's Mat.

Parameters
[in]aThe PETSc Mat object.
[in]refIf true, we increase the reference count of a.

Definition at line 826 of file petsc.cpp.

mfem::PetscParMatrix::PetscParMatrix ( const HypreParMatrix ha,
Operator::Type  tid 
)

Convert a HypreParMatrix ha to a PetscParMatrix in the given PETSc format tid.

The supported type ids are: Operator::PETSC_MATAIJ, Operator::PETSC_MATIS, and Operator::PETSC_MATSHELL.

Definition at line 371 of file petsc.cpp.

mfem::PetscParMatrix::PetscParMatrix ( MPI_Comm  comm,
const Operator op,
Operator::Type  tid 
)

Convert an mfem::Operator into a PetscParMatrix in the given PETSc format tid.

If tid is Operator::PETSC_MATSHELL and op is not a PetscParMatrix, it converts any mfem::Operator op implementing Operator::Mult() and Operator::MultTranspose() into a PetscParMatrix. The Operator op should not be deleted while the constructed PetscParMatrix is used.

Otherwise, it tries to convert the operator in PETSc's classes.

In particular, if op is a BlockOperator, then a MATNEST Mat object is created using tid as the type for the blocks.

Definition at line 390 of file petsc.cpp.

mfem::PetscParMatrix::PetscParMatrix ( MPI_Comm  comm,
PetscInt  glob_size,
PetscInt *  row_starts,
SparseMatrix diag,
Operator::Type  tid 
)

Creates block-diagonal square parallel matrix.

The block-diagonal is given by diag which must be in CSR format (finalized). The new PetscParMatrix does not take ownership of any of the input arrays. The type id tid can be either PETSC_MATAIJ (parallel distributed CSR) or PETSC_MATIS.

Definition at line 402 of file petsc.cpp.

mfem::PetscParMatrix::PetscParMatrix ( MPI_Comm  comm,
PetscInt  global_num_rows,
PetscInt  global_num_cols,
PetscInt *  row_starts,
PetscInt *  col_starts,
SparseMatrix diag,
Operator::Type  tid 
)

Creates block-diagonal rectangular parallel matrix.

The block-diagonal is given by diag which must be in CSR format (finalized). The new PetscParMatrix does not take ownership of any of the input arrays. The type id tid can be either PETSC_MATAIJ (parallel distributed CSR) or PETSC_MATIS.

Definition at line 414 of file petsc.cpp.

virtual mfem::PetscParMatrix::~PetscParMatrix ( )
inlinevirtual

Calls PETSc's destroy function.

Definition at line 210 of file petsc.hpp.

Member Function Documentation

void mfem::PetscParMatrix::ConvertOperator ( MPI_Comm  comm,
const Operator op,
Mat *  B,
bool  assembled 
)
protected

Convert an mfem::Operator into a Mat B; op can be destroyed.

Definition at line 635 of file petsc.cpp.

void mfem::PetscParMatrix::Destroy ( )
protected

Delete all owned data. Does not perform re-initialization with defaults.

Definition at line 813 of file petsc.cpp.

void mfem::PetscParMatrix::EliminateRowsCols ( const Array< int > &  rows_cols,
const HypreParVector X,
HypreParVector B 
)

Eliminate rows and columns from the matrix, and rows from the vector B. Modify B with the BC values in X.

Definition at line 1129 of file petsc.cpp.

PetscParMatrix * mfem::PetscParMatrix::EliminateRowsCols ( const Array< int > &  rows_cols)

Eliminate rows and columns from the matrix and store the eliminated elements in a new matrix Ae (returned).

The sum of the modified matrix and the returned matrix, Ae, is equal to the original matrix.

Definition at line 1107 of file petsc.cpp.

MPI_Comm mfem::PetscParMatrix::GetComm ( ) const
inline

Get the associated MPI communicator.

Definition at line 232 of file petsc.hpp.

PetscInt mfem::PetscParMatrix::GetGlobalNumCols ( ) const
inline

Returns the global number of columns.

Definition at line 253 of file petsc.hpp.

PetscInt mfem::PetscParMatrix::GetGlobalNumRows ( ) const
inline

Returns the global number of rows.

Definition at line 250 of file petsc.hpp.

PetscInt mfem::PetscParMatrix::GetNumCols ( ) const

Returns the local number of columns.

Definition at line 331 of file petsc.cpp.

PetscInt mfem::PetscParMatrix::GetNumRows ( ) const

Returns the local number of rows.

Definition at line 324 of file petsc.cpp.

Operator::Type mfem::PetscParMatrix::GetType ( ) const

Definition at line 1149 of file petsc.cpp.

PetscParVector * mfem::PetscParMatrix::GetX ( ) const

Returns the inner vector in the domain of A (it creates it if needed)

Definition at line 906 of file petsc.cpp.

PetscParVector * mfem::PetscParMatrix::GetY ( ) const

Returns the inner vector in the range of A (it creates it if needed)

Definition at line 916 of file petsc.cpp.

void mfem::PetscParMatrix::Init ( )
protected

Initialize with defaults. Does not initialize inherited members.

Definition at line 359 of file petsc.cpp.

PetscInt mfem::PetscParMatrix::M ( ) const

Returns the global number of rows.

Definition at line 338 of file petsc.cpp.

void mfem::PetscParMatrix::MakeRef ( const PetscParMatrix master)

Makes this object a reference to another PetscParMatrix.

Definition at line 896 of file petsc.cpp.

void mfem::PetscParMatrix::MakeWrapper ( MPI_Comm  comm,
const Operator op,
Mat *  B 
)
protected

Creates a wrapper around a mfem::Operator op using PETSc's MATSHELL object and returns the Mat in B.

This does not take any reference to op, that should not be destroyed until B is needed.

Definition at line 614 of file petsc.cpp.

void mfem::PetscParMatrix::Mult ( double  a,
const Vector x,
double  b,
Vector y 
) const

Matvec: y = a A x + b y.

Definition at line 945 of file petsc.cpp.

virtual void mfem::PetscParMatrix::Mult ( const Vector x,
Vector y 
) const
inlinevirtual

Operator application: y=A(x).

Implements mfem::Operator.

Definition at line 225 of file petsc.hpp.

void mfem::PetscParMatrix::MultTranspose ( double  a,
const Vector x,
double  b,
Vector y 
) const

Matvec transpose: y = a A^T x + b y.

Definition at line 961 of file petsc.cpp.

virtual void mfem::PetscParMatrix::MultTranspose ( const Vector x,
Vector y 
) const
inlinevirtual

Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an error.

Reimplemented from mfem::Operator.

Definition at line 228 of file petsc.hpp.

PetscInt mfem::PetscParMatrix::N ( ) const

Returns the global number of columns.

Definition at line 345 of file petsc.cpp.

PetscInt mfem::PetscParMatrix::NNZ ( ) const

Returns the number of nonzeros.

Differently from HYPRE, this call is collective on the communicator, as this number is not stored inside PETSc, but needs to be computed.

Definition at line 352 of file petsc.cpp.

mfem::PetscParMatrix::operator Mat ( ) const
inline

Typecasting to PETSc's Mat type.

Definition at line 235 of file petsc.hpp.

void mfem::PetscParMatrix::operator*= ( double  s)

Scale all entries by s: A_scaled = s*A.

Definition at line 940 of file petsc.cpp.

PetscParMatrix & mfem::PetscParMatrix::operator+= ( const PetscParMatrix B)

Definition at line 463 of file petsc.cpp.

PetscParMatrix & mfem::PetscParMatrix::operator= ( const PetscParMatrix B)

Definition at line 447 of file petsc.cpp.

PetscParMatrix & mfem::PetscParMatrix::operator= ( const HypreParMatrix B)

Definition at line 427 of file petsc.cpp.

void mfem::PetscParMatrix::Print ( const char *  fname = NULL,
bool  binary = false 
) const

Prints the matrix (to stdout if fname is NULL)

Definition at line 978 of file petsc.cpp.

Mat mfem::PetscParMatrix::ReleaseMat ( bool  dereference)

Release the PETSc Mat object. If dereference is true, decrement the refcount of the Mat object.

Definition at line 1136 of file petsc.cpp.

PetscParMatrix * mfem::PetscParMatrix::Transpose ( bool  action = false)

Returns the transpose of the PetscParMatrix.

If action is false, the new matrix is constructed with the PETSc function MatTranspose(). If action is true, then the matrix is not actually transposed. Instead, an object that behaves like the transpose is returned.

Definition at line 926 of file petsc.cpp.

Friends And Related Function Documentation

friend class PetscLinearSolver
friend

Definition at line 154 of file petsc.hpp.

friend class PetscPreconditioner
friend

Definition at line 155 of file petsc.hpp.

Member Data Documentation

Mat mfem::PetscParMatrix::A
protected

The actual PETSc object.

Definition at line 132 of file petsc.hpp.

PetscParVector* mfem::PetscParMatrix::X
mutableprotected

Auxiliary vectors for typecasting.

Definition at line 135 of file petsc.hpp.

PetscParVector * mfem::PetscParMatrix::Y
mutableprotected

Definition at line 135 of file petsc.hpp.


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