MFEM v4.7.0
Finite element discretization library
|
Wrapper for hypre's ParCSR matrix class. More...
#include <hypre.hpp>
Public Member Functions | |
HypreParMatrix () | |
An empty matrix to be used as a reference to an existing matrix. | |
void | WrapHypreParCSRMatrix (hypre_ParCSRMatrix *a, bool owner=true) |
Converts hypre's format to HypreParMatrix. | |
HypreParMatrix (hypre_ParCSRMatrix *a, bool owner=true) | |
Converts hypre's format to HypreParMatrix. | |
HypreParMatrix (MPI_Comm comm, HYPRE_BigInt glob_size, HYPRE_BigInt *row_starts, SparseMatrix *diag) | |
Creates block-diagonal square parallel matrix. | |
HypreParMatrix (MPI_Comm comm, HYPRE_BigInt global_num_rows, HYPRE_BigInt global_num_cols, HYPRE_BigInt *row_starts, HYPRE_BigInt *col_starts, SparseMatrix *diag) | |
Creates block-diagonal rectangular parallel matrix. | |
HypreParMatrix (MPI_Comm comm, HYPRE_BigInt global_num_rows, HYPRE_BigInt global_num_cols, HYPRE_BigInt *row_starts, HYPRE_BigInt *col_starts, SparseMatrix *diag, SparseMatrix *offd, HYPRE_BigInt *cmap, bool own_diag_offd=false) | |
Creates general (rectangular) parallel matrix. | |
HypreParMatrix (MPI_Comm comm, HYPRE_BigInt global_num_rows, HYPRE_BigInt global_num_cols, HYPRE_BigInt *row_starts, HYPRE_BigInt *col_starts, HYPRE_Int *diag_i, HYPRE_Int *diag_j, real_t *diag_data, HYPRE_Int *offd_i, HYPRE_Int *offd_j, real_t *offd_data, HYPRE_Int offd_num_cols, HYPRE_BigInt *offd_col_map, bool hypre_arrays=false) | |
Creates general (rectangular) parallel matrix. | |
HypreParMatrix (MPI_Comm comm, HYPRE_BigInt *row_starts, HYPRE_BigInt *col_starts, SparseMatrix *a) | |
Creates a parallel matrix from SparseMatrix on processor 0. | |
HypreParMatrix (MPI_Comm comm, HYPRE_BigInt global_num_rows, HYPRE_BigInt global_num_cols, HYPRE_BigInt *row_starts, HYPRE_BigInt *col_starts, Table *diag) | |
Creates boolean block-diagonal rectangular parallel matrix. | |
HypreParMatrix (MPI_Comm comm, int id, int np, HYPRE_BigInt *row, HYPRE_BigInt *col, HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd, HYPRE_Int *j_offd, HYPRE_BigInt *cmap, HYPRE_Int cmap_size) | |
Creates boolean rectangular parallel matrix. | |
HypreParMatrix (MPI_Comm comm, int nrows, HYPRE_BigInt glob_nrows, HYPRE_BigInt glob_ncols, int *I, HYPRE_BigInt *J, real_t *data, HYPRE_BigInt *rows, HYPRE_BigInt *cols) | |
Creates a general parallel matrix from a local CSR matrix on each processor described by the I, J and data arrays. | |
HypreParMatrix (const HypreParMatrix &P) | |
Copy constructor for a ParCSR matrix which creates a deep copy of structure and data from P. | |
void | MakeRef (const HypreParMatrix &master) |
Make this HypreParMatrix a reference to 'master'. | |
MPI_Comm | GetComm () const |
MPI communicator. | |
operator hypre_ParCSRMatrix * () const | |
Typecasting to hypre's hypre_ParCSRMatrix*. | |
operator HYPRE_ParCSRMatrix () | |
Typecasting to hypre's HYPRE_ParCSRMatrix, a.k.a. void *. | |
hypre_ParCSRMatrix * | StealData () |
Changes the ownership of the matrix. | |
void | SetOwnerFlags (signed char diag, signed char offd, signed char colmap) |
Explicitly set the three ownership flags, see docs for diagOwner etc. | |
signed char | OwnsDiag () const |
Get diag ownership flag. | |
signed char | OwnsOffd () const |
Get offd ownership flag. | |
signed char | OwnsColMap () const |
Get colmap ownership flag. | |
void | CopyRowStarts () |
void | CopyColStarts () |
HYPRE_BigInt | NNZ () const |
Returns the global number of nonzeros. | |
HYPRE_BigInt * | RowPart () |
Returns the row partitioning. | |
HYPRE_BigInt * | ColPart () |
Returns the column partitioning. | |
const HYPRE_BigInt * | RowPart () const |
Returns the row partitioning (const version) | |
const HYPRE_BigInt * | ColPart () const |
Returns the column partitioning (const version) | |
HYPRE_BigInt | M () const |
Returns the global number of rows. | |
HYPRE_BigInt | N () const |
Returns the global number of columns. | |
void | GetDiag (Vector &diag) const |
Get the local diagonal of the matrix. | |
void | GetDiag (SparseMatrix &diag) const |
Get the local diagonal block. NOTE: 'diag' will not own any data. | |
void | GetOffd (SparseMatrix &offd, HYPRE_BigInt *&cmap) const |
Get the local off-diagonal block. NOTE: 'offd' will not own any data. | |
void | MergeDiagAndOffd (SparseMatrix &merged) |
Get a single SparseMatrix containing all rows from this processor, merged from the diagonal and off-diagonal blocks stored by the HypreParMatrix. | |
void | AssembleDiagonal (Vector &diag) const override |
Return the diagonal of the matrix (Operator interface). | |
void | GetBlocks (Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const |
HypreParMatrix * | Transpose () const |
Returns the transpose of *this. | |
HypreParMatrix * | ExtractSubmatrix (const Array< int > &indices, real_t threshold=0.0) const |
int | GetNumRows () const |
Returns the number of rows in the diagonal block of the ParCSRMatrix. | |
int | GetNumCols () const |
Returns the number of columns in the diagonal block of the ParCSRMatrix. | |
HYPRE_BigInt | GetGlobalNumRows () const |
Return the global number of rows. | |
HYPRE_BigInt | GetGlobalNumCols () const |
Return the global number of columns. | |
HYPRE_BigInt * | GetRowStarts () const |
Return the parallel row partitioning array. | |
HYPRE_BigInt * | GetColStarts () const |
Return the parallel column partitioning array. | |
MemoryClass | GetMemoryClass () const override |
Return the MemoryClass preferred by the Operator. | |
void | EnsureMultTranspose () const |
Ensure the action of the transpose is performed fast. | |
void | ResetTranspose () const |
Reset (destroy) the internal transpose matrix that is created by EnsureMultTranspose() and MultTranspose(). | |
HYPRE_Int | Mult (HypreParVector &x, HypreParVector &y, real_t alpha=1.0, real_t beta=0.0) const |
Computes y = alpha * A * x + beta * y. | |
HYPRE_Int | Mult (HYPRE_ParVector x, HYPRE_ParVector y, real_t alpha=1.0, real_t beta=0.0) const |
Computes y = alpha * A * x + beta * y. | |
HYPRE_Int | MultTranspose (HypreParVector &x, HypreParVector &y, real_t alpha=1.0, real_t beta=0.0) const |
Computes y = alpha * A^t * x + beta * y. | |
void | Mult (real_t a, const Vector &x, real_t b, Vector &y) const |
void | MultTranspose (real_t a, const Vector &x, real_t b, Vector &y) const |
Computes y = alpha * A^t * x + beta * y. | |
void | Mult (const Vector &x, Vector &y) const override |
Operator application: y=A(x) . | |
void | MultTranspose (const Vector &x, Vector &y) const override |
Computes y = A^t * x. | |
void | AddMult (const Vector &x, Vector &y, const real_t a=1.0) const override |
Operator application: y+=A(x) (default) or y+=a*A(x) . | |
void | AddMultTranspose (const Vector &x, Vector &y, const real_t a=1.0) const override |
Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x) . | |
void | AbsMult (real_t a, const Vector &x, real_t b, Vector &y) const |
Computes y = a * |A| * x + b * y, using entry-wise absolute values of the matrix A. | |
void | AbsMultTranspose (real_t a, const Vector &x, real_t b, Vector &y) const |
Computes y = a * |At| * x + b * y, using entry-wise absolute values of the transpose of the matrix A. | |
void | BooleanMult (int alpha, const int *x, int beta, int *y) |
The "Boolean" analog of y = alpha * A * x + beta * y, where elements in the sparsity pattern of the matrix are treated as "true". | |
void | BooleanMultTranspose (int alpha, const int *x, int beta, int *y) |
The "Boolean" analog of y = alpha * A^T * x + beta * y, where elements in the sparsity pattern of the matrix are treated as "true". | |
HypreParMatrix & | operator= (real_t value) |
Initialize all entries with value. | |
HypreParMatrix & | operator+= (const HypreParMatrix &B) |
HypreParMatrix & | Add (const real_t beta, const HypreParMatrix &B) |
HypreParMatrix * | LeftDiagMult (const SparseMatrix &D, HYPRE_BigInt *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. | |
void | ScaleRows (const Vector &s) |
Scale the local row i by s(i). | |
void | InvScaleRows (const Vector &s) |
Scale the local row i by 1./s(i) | |
void | operator*= (real_t s) |
Scale all entries by s: A_scaled = s*A. | |
void | Threshold (real_t threshold=0.0) |
Remove values smaller in absolute value than some threshold. | |
void | DropSmallEntries (real_t tol) |
Wrapper for hypre_ParCSRMatrixDropSmallEntries in different versions of hypre. Drop off-diagonal entries that are smaller than tol * l2 norm of its row. | |
void | EliminateZeroRows () |
If a row contains only zeros, set its diagonal to 1. | |
void | EliminateRowsCols (const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B) |
HypreParMatrix * | EliminateRowsCols (const Array< int > &rows_cols) |
HypreParMatrix * | EliminateCols (const Array< int > &cols) |
void | EliminateRows (const Array< int > &rows) |
Eliminate rows from the diagonal and off-diagonal blocks of the matrix. | |
void | EliminateBC (const HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B) const |
Eliminate essential BC specified by ess_dof_list from the solution X to the r.h.s. B. | |
void | EliminateBC (const Array< int > &ess_dofs, DiagonalPolicy diag_policy) |
Eliminate essential (Dirichlet) boundary conditions. | |
void | HostRead () const |
Update the internal hypre_ParCSRMatrix object, A, to be on host. | |
void | HostReadWrite () |
Update the internal hypre_ParCSRMatrix object, A, to be on host. | |
void | HostWrite () |
Update the internal hypre_ParCSRMatrix object, A, to be on host. | |
void | HypreRead () const |
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space. | |
void | HypreReadWrite () |
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space. | |
void | HypreWrite () |
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space. | |
Memory< HYPRE_Int > & | GetDiagMemoryI () |
Memory< HYPRE_Int > & | GetDiagMemoryJ () |
Memory< real_t > & | GetDiagMemoryData () |
const Memory< HYPRE_Int > & | GetDiagMemoryI () const |
const Memory< HYPRE_Int > & | GetDiagMemoryJ () const |
const Memory< real_t > & | GetDiagMemoryData () const |
void | Print (const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0) const |
Prints the locally owned rows in parallel. | |
void | Read (MPI_Comm comm, const char *fname) |
Reads the matrix from a file. | |
void | Read_IJMatrix (MPI_Comm comm, const char *fname) |
Read a matrix saved as a HYPRE_IJMatrix. | |
void | PrintCommPkg (std::ostream &out=mfem::out) const |
Print information about the hypre_ParCSRCommPkg of the HypreParMatrix. | |
void | PrintHash (std::ostream &out) const |
Print sizes and hashes for all data arrays of the HypreParMatrix from the local MPI rank. | |
virtual | ~HypreParMatrix () |
Calls hypre's destroy function. | |
Type | GetType () const |
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 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 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. | |
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... | |
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::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. | |
mfem::HypreParMatrix::HypreParMatrix | ( | ) |
|
inlineexplicit |
Converts hypre's format to HypreParMatrix.
If owner is false, ownership of a is not transferred
mfem::HypreParMatrix::HypreParMatrix | ( | MPI_Comm | comm, |
HYPRE_BigInt | glob_size, | ||
HYPRE_BigInt * | 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. See here for a description of the row partitioning array row_starts.
mfem::HypreParMatrix::HypreParMatrix | ( | MPI_Comm | comm, |
HYPRE_BigInt | global_num_rows, | ||
HYPRE_BigInt | global_num_cols, | ||
HYPRE_BigInt * | row_starts, | ||
HYPRE_BigInt * | 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. See here for a description of the partitioning arrays row_starts and col_starts.
mfem::HypreParMatrix::HypreParMatrix | ( | MPI_Comm | comm, |
HYPRE_BigInt | global_num_rows, | ||
HYPRE_BigInt | global_num_cols, | ||
HYPRE_BigInt * | row_starts, | ||
HYPRE_BigInt * | col_starts, | ||
SparseMatrix * | diag, | ||
SparseMatrix * | offd, | ||
HYPRE_BigInt * | cmap, | ||
bool | own_diag_offd = false ) |
Creates general (rectangular) parallel matrix.
The new HypreParMatrix does not take ownership of any of the input arrays, if own_diag_offd is false (default). If own_diag_offd is true, ownership of diag and offd is transferred to the HypreParMatrix.
See here for a description of the partitioning arrays row_starts and col_starts.
mfem::HypreParMatrix::HypreParMatrix | ( | MPI_Comm | comm, |
HYPRE_BigInt | global_num_rows, | ||
HYPRE_BigInt | global_num_cols, | ||
HYPRE_BigInt * | row_starts, | ||
HYPRE_BigInt * | col_starts, | ||
HYPRE_Int * | diag_i, | ||
HYPRE_Int * | diag_j, | ||
real_t * | diag_data, | ||
HYPRE_Int * | offd_i, | ||
HYPRE_Int * | offd_j, | ||
real_t * | offd_data, | ||
HYPRE_Int | offd_num_cols, | ||
HYPRE_BigInt * | offd_col_map, | ||
bool | hypre_arrays = false ) |
Creates general (rectangular) parallel matrix.
The new HypreParMatrix takes ownership of all input arrays, except col_starts and row_starts. See here for a description of the partitioning arrays row_starts and col_starts.
If hypre_arrays is false, all arrays (except row_starts and col_starts) are assumed to be allocated according to the MemoryType returned by Device::GetHostMemoryType(). If hypre_arrays is true, then the same arrays are assumed to be allocated by hypre as host arrays.
mfem::HypreParMatrix::HypreParMatrix | ( | MPI_Comm | comm, |
HYPRE_BigInt * | row_starts, | ||
HYPRE_BigInt * | col_starts, | ||
SparseMatrix * | a ) |
Creates a parallel matrix from SparseMatrix on processor 0.
See here for a description of the partitioning arrays row_starts and col_starts.
mfem::HypreParMatrix::HypreParMatrix | ( | MPI_Comm | comm, |
HYPRE_BigInt | global_num_rows, | ||
HYPRE_BigInt | global_num_cols, | ||
HYPRE_BigInt * | row_starts, | ||
HYPRE_BigInt * | col_starts, | ||
Table * | diag ) |
Creates boolean block-diagonal rectangular parallel matrix.
The new HypreParMatrix does not take ownership of any of the input arrays. See here for a description of the partitioning arrays row_starts and col_starts.
mfem::HypreParMatrix::HypreParMatrix | ( | MPI_Comm | comm, |
int | id, | ||
int | np, | ||
HYPRE_BigInt * | row, | ||
HYPRE_BigInt * | col, | ||
HYPRE_Int * | i_diag, | ||
HYPRE_Int * | j_diag, | ||
HYPRE_Int * | i_offd, | ||
HYPRE_Int * | j_offd, | ||
HYPRE_BigInt * | 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. See here for a description of the partitioning arrays row and col.
mfem::HypreParMatrix::HypreParMatrix | ( | MPI_Comm | comm, |
int | nrows, | ||
HYPRE_BigInt | glob_nrows, | ||
HYPRE_BigInt | glob_ncols, | ||
int * | I, | ||
HYPRE_BigInt * | J, | ||
real_t * | data, | ||
HYPRE_BigInt * | rows, | ||
HYPRE_BigInt * | 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). See here for a description of the partitioning arrays rows and cols.
mfem::HypreParMatrix::HypreParMatrix | ( | const HypreParMatrix & | P | ) |
|
inlinevirtual |
|
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, or B has an empty off-diagonal block. We also assume that the sparsity pattern of *this
contains that of B
. For a more general case consider the stand-alone function ParAdd described below.
|
inlineoverridevirtual |
Operator application: y+=A(x)
(default) or y+=a*A(x)
.
Reimplemented from mfem::Operator.
|
inlineoverridevirtual |
Operator transpose application: y+=A^t(x)
(default) or y+=a*A^t(x)
.
Reimplemented from mfem::Operator.
|
inlineoverridevirtual |
Return the diagonal of the matrix (Operator interface).
Reimplemented from mfem::Operator.
|
inline |
|
inline |
|
inline |
|
inline |
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.
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.
void mfem::HypreParMatrix::DropSmallEntries | ( | real_t | tol | ) |
Wrapper for hypre_ParCSRMatrixDropSmallEntries in different versions of hypre. Drop off-diagonal entries that are smaller than tol * l2 norm of its row.
For HYPRE versions < 2.14, this method just calls Threshold() with threshold = tol * max(l2 row norm).
void mfem::HypreParMatrix::EliminateBC | ( | const Array< int > & | ess_dofs, |
DiagonalPolicy | diag_policy ) |
void mfem::HypreParMatrix::EliminateBC | ( | const HypreParMatrix & | Ae, |
const Array< int > & | ess_dof_list, | ||
const Vector & | X, | ||
Vector & | B ) const |
HypreParMatrix * mfem::HypreParMatrix::EliminateCols | ( | const Array< int > & | cols | ) |
void mfem::HypreParMatrix::EliminateRows | ( | const Array< int > & | rows | ) |
HypreParMatrix * mfem::HypreParMatrix::EliminateRowsCols | ( | const Array< int > & | rows_cols | ) |
void mfem::HypreParMatrix::EliminateRowsCols | ( | const Array< int > & | rows_cols, |
const HypreParVector & | X, | ||
HypreParVector & | B ) |
|
inline |
void mfem::HypreParMatrix::EnsureMultTranspose | ( | ) | const |
Ensure the action of the transpose is performed fast.
When HYPRE is built for GPUs, this method will construct and store the transposes of the 'diag' and 'offd' CSR matrices. When HYPRE is not built for GPUs, this method is a no-op.
This method is automatically called by MultTranspose().
If the matrix is modified the old transpose blocks can be deleted by calling ResetTranspose().
HypreParMatrix * mfem::HypreParMatrix::ExtractSubmatrix | ( | const Array< int > & | indices, |
real_t | threshold = 0.0 ) const |
void mfem::HypreParMatrix::GetBlocks | ( | Array2D< HypreParMatrix * > & | blocks, |
bool | interleaved_rows = false, | ||
bool | interleaved_cols = false ) const |
|
inline |
|
inline |
void mfem::HypreParMatrix::GetDiag | ( | SparseMatrix & | diag | ) | const |
void mfem::HypreParMatrix::GetDiag | ( | Vector & | diag | ) | const |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inlineoverridevirtual |
Return the MemoryClass preferred by the Operator.
This is the MemoryClass that will be used to access the input and output vectors in the Mult() and MultTranspose() methods.
For example, classes using the mfem::forall macro for implementation can return the value returned by Device::GetMemoryClass().
The default implementation of this method in class Operator returns MemoryClass::HOST.
Reimplemented from mfem::Operator.
|
inline |
|
inline |
void mfem::HypreParMatrix::GetOffd | ( | SparseMatrix & | offd, |
HYPRE_BigInt *& | cmap ) const |
|
inline |
|
inline |
|
inline |
Update the internal hypre_ParCSRMatrix object, A, to be on host.
After this call A's diagonal and off-diagonal can be modified on host and subsequent calls to Hypre{Read,Write,ReadWrite} will require a deep copy of the data if hypre is built with device support.
|
inline |
Update the internal hypre_ParCSRMatrix object, A, to be on host.
Similar to HostReadWrite(), except that the data will never be copied from device to host to ensure host contains the correct current data.
|
inline |
|
inline |
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
After this call A's diagonal and off-diagonal can be modified in hypre memory space and subsequent calls to Host{Read,Write,ReadWrite} will require a deep copy of the data if hypre is built with device support.
|
inline |
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
Similar to HostReadWrite(), except that the data will never be copied from host to hypre memory space to ensure the latter contains the correct current data.
void mfem::HypreParMatrix::InvScaleRows | ( | const Vector & | s | ) |
HypreParMatrix * mfem::HypreParMatrix::LeftDiagMult | ( | const SparseMatrix & | D, |
HYPRE_BigInt * | 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.
|
inline |
void mfem::HypreParMatrix::MakeRef | ( | const HypreParMatrix & | master | ) |
Make this HypreParMatrix a reference to 'master'.
void mfem::HypreParMatrix::MergeDiagAndOffd | ( | SparseMatrix & | merged | ) |
Get a single SparseMatrix containing all rows from this processor, merged from the diagonal and off-diagonal blocks stored by the HypreParMatrix.
HYPRE_Int mfem::HypreParMatrix::Mult | ( | HypreParVector & | x, |
HypreParVector & | y, | ||
real_t | alpha = 1.0, | ||
real_t | beta = 0.0 ) const |
|
inlineoverridevirtual |
Computes y = A^t * x.
If the matrix is modified, call ResetTranspose() and optionally EnsureMultTranspose() to make sure this method uses the correct updated transpose.
Reimplemented from mfem::Operator.
HYPRE_Int mfem::HypreParMatrix::MultTranspose | ( | HypreParVector & | x, |
HypreParVector & | y, | ||
real_t | alpha = 1.0, | ||
real_t | beta = 0.0 ) const |
Computes y = alpha * A^t * x + beta * y.
If the matrix is modified, call ResetTranspose() and optionally EnsureMultTranspose() to make sure this method uses the correct updated transpose.
Computes y = alpha * A^t * x + beta * y.
If the matrix is modified, call ResetTranspose() and optionally EnsureMultTranspose() to make sure this method uses the correct updated transpose.
|
inline |
|
inline |
|
inline |
|
inline |
void mfem::HypreParMatrix::operator*= | ( | real_t | s | ) |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
void mfem::HypreParMatrix::Print | ( | const char * | fname, |
HYPRE_Int | offi = 0, | ||
HYPRE_Int | offj = 0 ) const |
void mfem::HypreParMatrix::PrintCommPkg | ( | std::ostream & | out = mfem::out | ) | const |
Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
void mfem::HypreParMatrix::PrintHash | ( | std::ostream & | out | ) | const |
Print sizes and hashes for all data arrays of the HypreParMatrix from the local MPI rank.
This is a compact text representation of the local data of the HypreParMatrix that can be used to compare matrices from different runs without the need to save the whole matrix.
void mfem::HypreParMatrix::Read | ( | MPI_Comm | comm, |
const char * | fname ) |
void mfem::HypreParMatrix::Read_IJMatrix | ( | MPI_Comm | comm, |
const char * | fname ) |
void mfem::HypreParMatrix::ResetTranspose | ( | ) | const |
Reset (destroy) the internal transpose matrix that is created by EnsureMultTranspose() and MultTranspose().
If the matrix is modified, this method should be called to delete the out-of-date transpose that is stored internally.
|
inline |
|
inline |
void mfem::HypreParMatrix::ScaleRows | ( | const Vector & | s | ) |
void mfem::HypreParMatrix::SetOwnerFlags | ( | signed char | diag, |
signed char | offd, | ||
signed char | colmap ) |
hypre_ParCSRMatrix * mfem::HypreParMatrix::StealData | ( | ) |
void mfem::HypreParMatrix::Threshold | ( | real_t | threshold = 0.0 | ) |
HypreParMatrix * mfem::HypreParMatrix::Transpose | ( | ) | const |
void mfem::HypreParMatrix::WrapHypreParCSRMatrix | ( | hypre_ParCSRMatrix * | a, |
bool | owner = true ) |
Converts hypre's format to HypreParMatrix.
If owner is false, ownership of a is not transferred