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

Wrapper for hypre's ParCSR matrix class. More...

#include <hypre.hpp>

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

Public Member Functions

 HypreParMatrix ()
 An empty matrix to be used as a reference to an existing matrix. More...
 
 HypreParMatrix (hypre_ParCSRMatrix *a, bool owner=true)
 Converts hypre's format to HypreParMatrix. More...
 
 HypreParMatrix (MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *row_starts, SparseMatrix *diag)
 
 HypreParMatrix (MPI_Comm comm, HYPRE_Int global_num_rows, HYPRE_Int global_num_cols, HYPRE_Int *row_starts, HYPRE_Int *col_starts, SparseMatrix *diag)
 
 HypreParMatrix (MPI_Comm comm, HYPRE_Int global_num_rows, HYPRE_Int global_num_cols, HYPRE_Int *row_starts, HYPRE_Int *col_starts, SparseMatrix *diag, SparseMatrix *offd, HYPRE_Int *cmap)
 
 HypreParMatrix (MPI_Comm comm, HYPRE_Int global_num_rows, HYPRE_Int global_num_cols, HYPRE_Int *row_starts, HYPRE_Int *col_starts, HYPRE_Int *diag_i, HYPRE_Int *diag_j, double *diag_data, HYPRE_Int *offd_i, HYPRE_Int *offd_j, double *offd_data, HYPRE_Int offd_num_cols, HYPRE_Int *offd_col_map)
 
 HypreParMatrix (MPI_Comm comm, HYPRE_Int *row_starts, HYPRE_Int *col_starts, SparseMatrix *a)
 Creates a parallel matrix from SparseMatrix on processor 0. More...
 
 HypreParMatrix (MPI_Comm comm, HYPRE_Int global_num_rows, HYPRE_Int global_num_cols, HYPRE_Int *row_starts, HYPRE_Int *col_starts, Table *diag)
 
 HypreParMatrix (MPI_Comm comm, int id, int np, HYPRE_Int *row, HYPRE_Int *col, HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd, HYPRE_Int *j_offd, HYPRE_Int *cmap, HYPRE_Int cmap_size)
 
 HypreParMatrix (MPI_Comm comm, int nrows, HYPRE_Int glob_nrows, HYPRE_Int glob_ncols, int *I, HYPRE_Int *J, double *data, HYPRE_Int *rows, HYPRE_Int *cols)
 
void MakeRef (const HypreParMatrix &master)
 Make this HypreParMatrix a reference to 'master'. More...
 
MPI_Comm GetComm () const
 MPI communicator. More...
 
 operator hypre_ParCSRMatrix * () const
 Typecasting to hypre's hypre_ParCSRMatrix*. More...
 
 operator HYPRE_ParCSRMatrix ()
 Typecasting to hypre's HYPRE_ParCSRMatrix, a.k.a. void *. More...
 
hypre_ParCSRMatrix * StealData ()
 Changes the ownership of the the matrix. More...
 
void SetOwnerFlags (signed char diag, signed char offd, signed char colmap)
 Explicitly set the three ownership flags, see docs for diagOwner etc. More...
 
signed char OwnsDiag () const
 Get diag ownership flag. More...
 
signed char OwnsOffd () const
 Get offd ownership flag. More...
 
signed char OwnsColMap () const
 Get colmap ownership flag. More...
 
void CopyRowStarts ()
 
void CopyColStarts ()
 
HYPRE_Int NNZ () const
 Returns the global number of nonzeros. More...
 
HYPRE_Int * RowPart ()
 Returns the row partitioning. More...
 
HYPRE_Int * ColPart ()
 Returns the column partitioning. More...
 
const HYPRE_Int * RowPart () const
 Returns the row partitioning (const version) More...
 
const HYPRE_Int * ColPart () const
 Returns the column partitioning (const version) More...
 
HYPRE_Int M () const
 Returns the global number of rows. More...
 
HYPRE_Int N () const
 Returns the global number of columns. More...
 
void GetDiag (Vector &diag) const
 Get the local diagonal of the matrix. More...
 
void GetDiag (SparseMatrix &diag) const
 Get the local diagonal block. NOTE: 'diag' will not own any data. More...
 
void GetOffd (SparseMatrix &offd, HYPRE_Int *&cmap) const
 Get the local off-diagonal block. NOTE: 'offd' will not own any data. More...
 
void GetBlocks (Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
 
HypreParMatrixTranspose () const
 Returns the transpose of *this. More...
 
int GetNumRows () const
 Returns the number of rows in the diagonal block of the ParCSRMatrix. More...
 
int GetNumCols () const
 Returns the number of columns in the diagonal block of the ParCSRMatrix. More...
 
HYPRE_Int GetGlobalNumRows () const
 
HYPRE_Int GetGlobalNumCols () const
 
HYPRE_Int * GetRowStarts () const
 
HYPRE_Int * GetColStarts () const
 
HYPRE_Int Mult (HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
 Computes y = alpha * A * x + beta * y. More...
 
HYPRE_Int Mult (HYPRE_ParVector x, HYPRE_ParVector y, double alpha=1.0, double beta=0.0)
 Computes y = alpha * A * x + beta * y. More...
 
HYPRE_Int MultTranspose (HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
 Computes y = alpha * A^t * x + beta * y. More...
 
void Mult (double a, const Vector &x, double b, Vector &y) const
 
void MultTranspose (double a, const Vector &x, double b, Vector &y) const
 
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...
 
void BooleanMult (int alpha, int *x, int beta, int *y)
 
void BooleanMultTranspose (int alpha, int *x, int beta, int *y)
 
HypreParMatrixoperator= (double value)
 Initialize all entries with value. More...
 
HypreParMatrixoperator+= (const HypreParMatrix &B)
 
HypreParMatrixAdd (const double beta, const HypreParMatrix &B)
 
HypreParMatrixLeftDiagMult (const SparseMatrix &D, HYPRE_Int *row_starts=NULL) const
 Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result as a new HypreParMatrix. More...
 
void ScaleRows (const Vector &s)
 Scale the local row i by s(i). More...
 
void InvScaleRows (const Vector &s)
 Scale the local row i by 1./s(i) More...
 
void operator*= (double s)
 Scale all entries by s: A_scaled = s*A. More...
 
void Threshold (double threshold=0.0)
 Remove values smaller in absolute value than some threshold. More...
 
void EliminateZeroRows ()
 If a row contains only zeros, set its diagonal to 1. More...
 
void EliminateRowsCols (const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
 
HypreParMatrixEliminateRowsCols (const Array< int > &rows_cols)
 
void Print (const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0)
 Prints the locally owned rows in parallel. More...
 
void Read (MPI_Comm comm, const char *fname)
 Reads the matrix from a file. More...
 
void Read_IJMatrix (MPI_Comm comm, const char *fname)
 Read a matrix saved as a HYPRE_IJMatrix. More...
 
void PrintCommPkg (std::ostream &out=mfem::out) const
 Print information about the hypre_ParCSRCommPkg of the HypreParMatrix. More...
 
virtual ~HypreParMatrix ()
 Calls hypre's destroy function. More...
 
Type GetType () const
 
- 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...
 
Type GetType () const
 Return the type ID of the Operator class. More...
 

Additional Inherited Members

- Public Types inherited from mfem::Operator
enum  Type {
  ANY_TYPE, MFEM_SPARSEMAT, Hypre_ParCSR, PETSC_MATAIJ,
  PETSC_MATIS, PETSC_MATSHELL, PETSC_MATNEST, PETSC_MATHYPRE,
  PETSC_MATGENERIC
}
 Enumeration defining IDs for some classes derived from Operator. More...
 
- 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...
 

Detailed Description

Wrapper for hypre's ParCSR matrix class.

Definition at line 175 of file hypre.hpp.

Constructor & Destructor Documentation

mfem::HypreParMatrix::HypreParMatrix ( )

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

Definition at line 304 of file hypre.cpp.

mfem::HypreParMatrix::HypreParMatrix ( hypre_ParCSRMatrix *  a,
bool  owner = true 
)
inlineexplicit

Converts hypre's format to HypreParMatrix.

If owner is false, ownership of a is not transferred

Definition at line 228 of file hypre.hpp.

mfem::HypreParMatrix::HypreParMatrix ( MPI_Comm  comm,
HYPRE_Int  glob_size,
HYPRE_Int *  row_starts,
SparseMatrix diag 
)

Creates block-diagonal square parallel matrix. Diagonal is given by diag which must be in CSR format (finalized). The new HypreParMatrix does not take ownership of any of the input arrays.

Warning
The ordering of the columns in each row in *diag may be changed by this constructor to ensure that the first entry in each row is the diagonal one. This is expected by most hypre functions.

Definition at line 362 of file hypre.cpp.

mfem::HypreParMatrix::HypreParMatrix ( MPI_Comm  comm,
HYPRE_Int  global_num_rows,
HYPRE_Int  global_num_cols,
HYPRE_Int *  row_starts,
HYPRE_Int *  col_starts,
SparseMatrix diag 
)

Creates block-diagonal rectangular parallel matrix. Diagonal is given by diag which must be in CSR format (finalized). The new HypreParMatrix does not take ownership of any of the input arrays.

Definition at line 397 of file hypre.cpp.

mfem::HypreParMatrix::HypreParMatrix ( MPI_Comm  comm,
HYPRE_Int  global_num_rows,
HYPRE_Int  global_num_cols,
HYPRE_Int *  row_starts,
HYPRE_Int *  col_starts,
SparseMatrix diag,
SparseMatrix offd,
HYPRE_Int *  cmap 
)

Creates general (rectangular) parallel matrix. The new HypreParMatrix does not take ownership of any of the input arrays.

Definition at line 434 of file hypre.cpp.

mfem::HypreParMatrix::HypreParMatrix ( MPI_Comm  comm,
HYPRE_Int  global_num_rows,
HYPRE_Int  global_num_cols,
HYPRE_Int *  row_starts,
HYPRE_Int *  col_starts,
HYPRE_Int *  diag_i,
HYPRE_Int *  diag_j,
double *  diag_data,
HYPRE_Int *  offd_i,
HYPRE_Int *  offd_j,
double *  offd_data,
HYPRE_Int  offd_num_cols,
HYPRE_Int *  offd_col_map 
)

Creates general (rectangular) parallel matrix. The new HypreParMatrix takes ownership of all input arrays, except col_starts and row_starts.

Definition at line 478 of file hypre.cpp.

mfem::HypreParMatrix::HypreParMatrix ( MPI_Comm  comm,
HYPRE_Int *  row_starts,
HYPRE_Int *  col_starts,
SparseMatrix a 
)

Creates a parallel matrix from SparseMatrix on processor 0.

Definition at line 532 of file hypre.cpp.

mfem::HypreParMatrix::HypreParMatrix ( MPI_Comm  comm,
HYPRE_Int  global_num_rows,
HYPRE_Int  global_num_cols,
HYPRE_Int *  row_starts,
HYPRE_Int *  col_starts,
Table diag 
)

Creates boolean block-diagonal rectangular parallel matrix. The new HypreParMatrix does not take ownership of any of the input arrays.

Definition at line 572 of file hypre.cpp.

mfem::HypreParMatrix::HypreParMatrix ( MPI_Comm  comm,
int  id,
int  np,
HYPRE_Int *  row,
HYPRE_Int *  col,
HYPRE_Int *  i_diag,
HYPRE_Int *  j_diag,
HYPRE_Int *  i_offd,
HYPRE_Int *  j_offd,
HYPRE_Int *  cmap,
HYPRE_Int  cmap_size 
)

Creates boolean rectangular parallel matrix. The new HypreParMatrix takes ownership of the arrays i_diag, j_diag, i_offd, j_offd, and cmap; does not take ownership of the arrays row and col.

Definition at line 612 of file hypre.cpp.

mfem::HypreParMatrix::HypreParMatrix ( MPI_Comm  comm,
int  nrows,
HYPRE_Int  glob_nrows,
HYPRE_Int  glob_ncols,
int *  I,
HYPRE_Int *  J,
double *  data,
HYPRE_Int *  rows,
HYPRE_Int *  cols 
)

Creates a general parallel matrix from a local CSR matrix on each processor described by the I, J and data arrays. The local matrix should be of size (local) nrows by (global) glob_ncols. The new parallel matrix contains copies of all input arrays (so they can be deleted).

Definition at line 693 of file hypre.cpp.

virtual mfem::HypreParMatrix::~HypreParMatrix ( )
inlinevirtual

Calls hypre's destroy function.

Definition at line 485 of file hypre.hpp.

Member Function Documentation

HypreParMatrix& mfem::HypreParMatrix::Add ( const double  beta,
const HypreParMatrix B 
)
inline

Perform the operation *this += beta*B, assuming that both matrices use the same row and column partitions and the same col_map_offd arrays. We also assume that the sparsity pattern of *this contains that of B.

Definition at line 433 of file hypre.hpp.

void mfem::HypreParMatrix::BooleanMult ( int  alpha,
int *  x,
int  beta,
int *  y 
)
inline

The "Boolean" analog of y = alpha * A * x + beta * y, where elements in the sparsity pattern of the matrix are treated as "true".

Definition at line 409 of file hypre.hpp.

void mfem::HypreParMatrix::BooleanMultTranspose ( int  alpha,
int *  x,
int  beta,
int *  y 
)
inline

The "Boolean" analog of y = alpha * A^T * x + beta * y, where elements in the sparsity pattern of the matrix are treated as "true".

Definition at line 416 of file hypre.hpp.

HYPRE_Int* mfem::HypreParMatrix::ColPart ( )
inline

Returns the column partitioning.

Definition at line 339 of file hypre.hpp.

const HYPRE_Int* mfem::HypreParMatrix::ColPart ( ) const
inline

Returns the column partitioning (const version)

Definition at line 343 of file hypre.hpp.

void mfem::HypreParMatrix::CopyColStarts ( )

If the HypreParMatrix does not own the col-starts array, make a copy of it that the HypreParMatrix will own. If the row-starts array is the same as the col-starts array, row-starts is also replaced.

Definition at line 882 of file hypre.cpp.

void mfem::HypreParMatrix::CopyRowStarts ( )

If the HypreParMatrix does not own the row-starts array, make a copy of it that the HypreParMatrix will own. If the col-starts array is the same as the row-starts array, col-starts is also replaced.

Definition at line 845 of file hypre.cpp.

void mfem::HypreParMatrix::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 1368 of file hypre.cpp.

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

Eliminate rows and columns from the matrix and store the eliminated elements in a new matrix Ae (returned), so that the modified matrix and Ae sum to the original matrix.

Definition at line 1379 of file hypre.cpp.

void mfem::HypreParMatrix::EliminateZeroRows ( )
inline

If a row contains only zeros, set its diagonal to 1.

Definition at line 462 of file hypre.hpp.

void mfem::HypreParMatrix::GetBlocks ( Array2D< HypreParMatrix * > &  blocks,
bool  interleaved_rows = false,
bool  interleaved_cols = false 
) const

Split the matrix into M x N equally sized blocks of parallel matrices. The size of 'blocks' must already be set to M x N.

Definition at line 965 of file hypre.cpp.

HYPRE_Int* mfem::HypreParMatrix::GetColStarts ( ) const
inline

Definition at line 387 of file hypre.hpp.

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

MPI communicator.

Definition at line 303 of file hypre.hpp.

void mfem::HypreParMatrix::GetDiag ( Vector diag) const

Get the local diagonal of the matrix.

Definition at line 923 of file hypre.cpp.

void mfem::HypreParMatrix::GetDiag ( SparseMatrix diag) const

Get the local diagonal block. NOTE: 'diag' will not own any data.

Definition at line 954 of file hypre.cpp.

HYPRE_Int mfem::HypreParMatrix::GetGlobalNumCols ( ) const
inline

Definition at line 382 of file hypre.hpp.

HYPRE_Int mfem::HypreParMatrix::GetGlobalNumRows ( ) const
inline

Definition at line 379 of file hypre.hpp.

int mfem::HypreParMatrix::GetNumCols ( ) const
inline

Returns the number of columns in the diagonal block of the ParCSRMatrix.

Definition at line 373 of file hypre.hpp.

int mfem::HypreParMatrix::GetNumRows ( ) const
inline

Returns the number of rows in the diagonal block of the ParCSRMatrix.

Definition at line 366 of file hypre.hpp.

void mfem::HypreParMatrix::GetOffd ( SparseMatrix offd,
HYPRE_Int *&  cmap 
) const

Get the local off-diagonal block. NOTE: 'offd' will not own any data.

Definition at line 959 of file hypre.cpp.

HYPRE_Int* mfem::HypreParMatrix::GetRowStarts ( ) const
inline

Definition at line 385 of file hypre.hpp.

Type mfem::HypreParMatrix::GetType ( ) const
inline

Definition at line 487 of file hypre.hpp.

void mfem::HypreParMatrix::InvScaleRows ( const Vector s)

Scale the local row i by 1./s(i)

Definition at line 1213 of file hypre.cpp.

HypreParMatrix * mfem::HypreParMatrix::LeftDiagMult ( const SparseMatrix D,
HYPRE_Int *  row_starts = NULL 
) const

Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result as a new HypreParMatrix.

If D has a different number of rows than A (this matrix), D's row starts array needs to be given (as returned by the methods GetDofOffsets/GetTrueDofOffsets of ParFiniteElementSpace). The new matrix D*A uses copies of the row- and column-starts arrays, so "this" matrix and row_starts can be deleted.

Note
This operation is local and does not require communication.

Definition at line 1081 of file hypre.cpp.

HYPRE_Int mfem::HypreParMatrix::M ( ) const
inline

Returns the global number of rows.

Definition at line 345 of file hypre.hpp.

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

Make this HypreParMatrix a reference to 'master'.

Definition at line 821 of file hypre.cpp.

HYPRE_Int mfem::HypreParMatrix::Mult ( HypreParVector x,
HypreParVector y,
double  alpha = 1.0,
double  beta = 0.0 
)

Computes y = alpha * A * x + beta * y.

Definition at line 1005 of file hypre.cpp.

HYPRE_Int mfem::HypreParMatrix::Mult ( HYPRE_ParVector  x,
HYPRE_ParVector  y,
double  alpha = 1.0,
double  beta = 0.0 
)

Computes y = alpha * A * x + beta * y.

Definition at line 1068 of file hypre.cpp.

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

Definition at line 1011 of file hypre.cpp.

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

Operator application: y=A(x).

Implements mfem::Operator.

Definition at line 402 of file hypre.hpp.

HYPRE_Int mfem::HypreParMatrix::MultTranspose ( HypreParVector x,
HypreParVector y,
double  alpha = 1.0,
double  beta = 0.0 
)

Computes y = alpha * A^t * x + beta * y.

Definition at line 1075 of file hypre.cpp.

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

Definition at line 1038 of file hypre.cpp.

virtual void mfem::HypreParMatrix::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 404 of file hypre.hpp.

HYPRE_Int mfem::HypreParMatrix::N ( ) const
inline

Returns the global number of columns.

Definition at line 347 of file hypre.hpp.

HYPRE_Int mfem::HypreParMatrix::NNZ ( ) const
inline

Returns the global number of nonzeros.

Definition at line 335 of file hypre.hpp.

mfem::HypreParMatrix::operator HYPRE_ParCSRMatrix ( )
inline

Typecasting to hypre's HYPRE_ParCSRMatrix, a.k.a. void *.

Definition at line 309 of file hypre.hpp.

mfem::HypreParMatrix::operator hypre_ParCSRMatrix * ( ) const
inline

Typecasting to hypre's hypre_ParCSRMatrix*.

Definition at line 306 of file hypre.hpp.

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

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

Definition at line 1254 of file hypre.cpp.

HypreParMatrix& mfem::HypreParMatrix::operator+= ( const HypreParMatrix B)
inline

Perform the operation *this += B, assuming that both matrices use the same row and column partitions and the same col_map_offd arrays. We also assume that the sparsity pattern of *this contains that of B.

Definition at line 428 of file hypre.hpp.

HypreParMatrix& mfem::HypreParMatrix::operator= ( double  value)
inline

Initialize all entries with value.

Definition at line 422 of file hypre.hpp.

signed char mfem::HypreParMatrix::OwnsColMap ( ) const
inline

Get colmap ownership flag.

Definition at line 323 of file hypre.hpp.

signed char mfem::HypreParMatrix::OwnsDiag ( ) const
inline

Get diag ownership flag.

Definition at line 319 of file hypre.hpp.

signed char mfem::HypreParMatrix::OwnsOffd ( ) const
inline

Get offd ownership flag.

Definition at line 321 of file hypre.hpp.

void mfem::HypreParMatrix::Print ( const char *  fname,
HYPRE_Int  offi = 0,
HYPRE_Int  offj = 0 
)

Prints the locally owned rows in parallel.

Definition at line 1391 of file hypre.cpp.

void mfem::HypreParMatrix::PrintCommPkg ( std::ostream &  out = mfem::out) const

Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.

Definition at line 1432 of file hypre.cpp.

void mfem::HypreParMatrix::Read ( MPI_Comm  comm,
const char *  fname 
)

Reads the matrix from a file.

Definition at line 1396 of file hypre.cpp.

void mfem::HypreParMatrix::Read_IJMatrix ( MPI_Comm  comm,
const char *  fname 
)

Read a matrix saved as a HYPRE_IJMatrix.

Definition at line 1411 of file hypre.cpp.

HYPRE_Int* mfem::HypreParMatrix::RowPart ( )
inline

Returns the row partitioning.

Definition at line 337 of file hypre.hpp.

const HYPRE_Int* mfem::HypreParMatrix::RowPart ( ) const
inline

Returns the row partitioning (const version)

Definition at line 341 of file hypre.hpp.

void mfem::HypreParMatrix::ScaleRows ( const Vector s)

Scale the local row i by s(i).

Definition at line 1178 of file hypre.cpp.

void mfem::HypreParMatrix::SetOwnerFlags ( signed char  diag,
signed char  offd,
signed char  colmap 
)
inline

Explicitly set the three ownership flags, see docs for diagOwner etc.

Definition at line 315 of file hypre.hpp.

hypre_ParCSRMatrix * mfem::HypreParMatrix::StealData ( )

Changes the ownership of the the matrix.

Definition at line 831 of file hypre.cpp.

void mfem::HypreParMatrix::Threshold ( double  threshold = 0.0)

Remove values smaller in absolute value than some threshold.

Definition at line 1292 of file hypre.cpp.

HypreParMatrix * mfem::HypreParMatrix::Transpose ( ) const

Returns the transpose of *this.

Definition at line 987 of file hypre.cpp.


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