MFEM v4.7.0
Finite element discretization library
|
Data type sparse matrix. More...
#include <sparsemat.hpp>
Public Member Functions | |
SparseMatrix () | |
Create an empty SparseMatrix. | |
SparseMatrix (int nrows, int ncols=-1) | |
Create a sparse matrix with flexible sparsity structure using a row-wise linked list (LIL) format. | |
SparseMatrix (int *i, int *j, real_t *data, int m, int n) | |
Create a sparse matrix in CSR format. Ownership of i, j, and data is transferred to the SparseMatrix. | |
SparseMatrix (int *i, int *j, real_t *data, int m, int n, bool ownij, bool owna, bool issorted) | |
Create a sparse matrix in CSR format. Ownership of i, j, and data is optionally transferred to the SparseMatrix. | |
SparseMatrix (int nrows, int ncols, int rowsize) | |
Create a sparse matrix in CSR format where each row has space allocated for exactly rowsize entries. | |
SparseMatrix (const SparseMatrix &mat, bool copy_graph=true, MemoryType mt=MemoryType::PRESERVE) | |
Copy constructor (deep copy). | |
SparseMatrix (const Vector &v) | |
Create a SparseMatrix with diagonal v, i.e. A = Diag(v) | |
void | OverrideSize (int height_, int width_) |
Sets the height and width of the matrix. | |
void | UseGPUSparse (bool useGPUSparse_=true) |
Runtime option to use cuSPARSE or hipSPARSE. Only valid when using a CUDA or HIP backend. | |
MFEM_DEPRECATED void | UseCuSparse (bool useCuSparse_=true) |
Deprecated equivalent of UseGPUSparse(). | |
SparseMatrix & | operator= (const SparseMatrix &rhs) |
Assignment operator: deep copy. | |
void | MakeRef (const SparseMatrix &master) |
Clear the contents of the SparseMatrix and make it a reference to master. | |
int | Size () const |
For backward compatibility, define Size() to be synonym of Height(). | |
void | Clear () |
Clear the contents of the SparseMatrix. | |
void | ClearGPUSparse () |
Clear the cuSPARSE/hipSPARSE descriptors. This must be called after releasing the device memory of A. | |
MFEM_DEPRECATED void | ClearCuSparse () |
Deprecated equivalent of ClearGPUSparse(). | |
bool | Empty () const |
Check if the SparseMatrix is empty. | |
int * | GetI () |
Return the array I. | |
const int * | GetI () const |
Return the array I, const version. | |
int * | GetJ () |
Return the array J. | |
const int * | GetJ () const |
Return the array J, const version. | |
real_t * | GetData () |
Return the element data, i.e. the array A. | |
const real_t * | GetData () const |
Return the element data, i.e. the array A, const version. | |
Memory< int > & | GetMemoryI () |
const Memory< int > & | GetMemoryI () const |
const int * | ReadI (bool on_dev=true) const |
int * | WriteI (bool on_dev=true) |
int * | ReadWriteI (bool on_dev=true) |
const int * | HostReadI () const |
int * | HostWriteI () |
int * | HostReadWriteI () |
Memory< int > & | GetMemoryJ () |
const Memory< int > & | GetMemoryJ () const |
const int * | ReadJ (bool on_dev=true) const |
int * | WriteJ (bool on_dev=true) |
int * | ReadWriteJ (bool on_dev=true) |
const int * | HostReadJ () const |
int * | HostWriteJ () |
int * | HostReadWriteJ () |
Memory< real_t > & | GetMemoryData () |
const Memory< real_t > & | GetMemoryData () const |
const real_t * | ReadData (bool on_dev=true) const |
real_t * | WriteData (bool on_dev=true) |
real_t * | ReadWriteData (bool on_dev=true) |
const real_t * | HostReadData () const |
real_t * | HostWriteData () |
real_t * | HostReadWriteData () |
int | RowSize (const int i) const |
Returns the number of elements in row i. | |
int | MaxRowSize () const |
Returns the maximum number of elements among all rows. | |
int * | GetRowColumns (const int row) |
Return a pointer to the column indices in a row. | |
const int * | GetRowColumns (const int row) const |
Return a pointer to the column indices in a row, const version. | |
real_t * | GetRowEntries (const int row) |
Return a pointer to the entries in a row. | |
const real_t * | GetRowEntries (const int row) const |
Return a pointer to the entries in a row, const version. | |
void | SetWidth (int width_=-1) |
Change the width of a SparseMatrix. | |
int | ActualWidth () const |
Returns the actual Width of the matrix. | |
void | SortColumnIndices () |
Sort the column indices corresponding to each row. | |
void | MoveDiagonalFirst () |
Move the diagonal entry to the first position in each row, preserving the order of the rest of the columns. | |
virtual real_t & | Elem (int i, int j) |
Returns reference to a_{ij}. | |
virtual const real_t & | Elem (int i, int j) const |
Returns constant reference to a_{ij}. | |
real_t & | operator() (int i, int j) |
Returns reference to A[i][j]. | |
const real_t & | operator() (int i, int j) const |
Returns reference to A[i][j]. | |
void | GetDiag (Vector &d) const |
Returns the Diagonal of A. | |
DenseMatrix * | ToDenseMatrix () const |
Produces a DenseMatrix from a SparseMatrix. | |
void | ToDenseMatrix (DenseMatrix &B) const |
Produces a DenseMatrix from a SparseMatrix. | |
virtual MemoryClass | GetMemoryClass () const |
Return the MemoryClass preferred by the Operator. | |
virtual void | Mult (const Vector &x, Vector &y) const |
Matrix vector multiplication. | |
virtual void | AddMult (const Vector &x, Vector &y, const real_t a=1.0) const |
y += A * x (default) or y += a * A * x | |
virtual void | MultTranspose (const Vector &x, Vector &y) const |
Multiply a vector with the transposed matrix. y = At * x. | |
virtual void | AddMultTranspose (const Vector &x, Vector &y, const real_t a=1.0) const |
y += At * x (default) or y += a * At * x | |
void | BuildTranspose () const |
Build and store internally the transpose of this matrix which will be used in the methods AddMultTranspose(), MultTranspose(), and AbsMultTranspose(). | |
void | ResetTranspose () const |
void | EnsureMultTranspose () const |
Ensures that the matrix is capable of performing MultTranspose(), AddMultTranspose(), and AbsMultTranspose(). | |
void | PartMult (const Array< int > &rows, const Vector &x, Vector &y) const |
void | PartAddMult (const Array< int > &rows, const Vector &x, Vector &y, const real_t a=1.0) const |
void | BooleanMult (const Array< int > &x, Array< int > &y) const |
y = A * x, treating all entries as booleans (zero=false, nonzero=true). | |
void | BooleanMultTranspose (const Array< int > &x, Array< int > &y) const |
y = At * x, treating all entries as booleans (zero=false, nonzero=true). | |
void | AbsMult (const Vector &x, Vector &y) const |
y = |A| * x, using entry-wise absolute values of matrix A | |
void | AbsMultTranspose (const Vector &x, Vector &y) const |
y = |At| * x, using entry-wise absolute values of the transpose of matrix A | |
real_t | InnerProduct (const Vector &x, const Vector &y) const |
Compute y^t A x. | |
void | GetRowSums (Vector &x) const |
For all i compute \( x_i = \sum_j A_{ij} \). | |
real_t | GetRowNorml1 (int irow) const |
For i = irow compute \( x_i = \sum_j | A_{i, j} | \). | |
virtual MatrixInverse * | Inverse () const |
This virtual method is not supported: it always returns NULL. | |
void | EliminateRow (int row, const real_t sol, Vector &rhs) |
Eliminates a column from the transpose matrix. | |
void | EliminateRow (int row, DiagonalPolicy dpolicy=DIAG_ZERO) |
Eliminates a row from the matrix. | |
void | EliminateCol (int col, DiagonalPolicy dpolicy=DIAG_ZERO) |
Eliminates the column col from the matrix. | |
void | EliminateCols (const Array< int > &cols, const Vector *x=NULL, Vector *b=NULL) |
Eliminate all columns i for which cols[i] != 0. | |
void | EliminateCols (const Array< int > &col_marker, SparseMatrix &Ae) |
Similar to EliminateCols + save the eliminated entries into Ae so that (*this) + Ae is equal to the original matrix. | |
void | EliminateRowCol (int rc, const real_t sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE) |
Eliminate row rc and column rc and modify the rhs using sol. | |
void | EliminateRowColMultipleRHS (int rc, const Vector &sol, DenseMatrix &rhs, DiagonalPolicy dpolicy=DIAG_ONE) |
Similar to EliminateRowCol(int, const double, Vector &, DiagonalPolicy), but multiple values for eliminated unknowns are accepted, and accordingly multiple right-hand-sides are used. | |
void | EliminateRowColDiag (int rc, real_t value) |
Perform elimination and set the diagonal entry to the given value. | |
void | EliminateRowCol (int rc, DiagonalPolicy dpolicy=DIAG_ONE) |
Eliminate row rc and column rc. | |
void | EliminateRowCol (int rc, SparseMatrix &Ae, DiagonalPolicy dpolicy=DIAG_ONE) |
Similar to EliminateRowCol(int, DiagonalPolicy) + save the eliminated entries into Ae so that (*this) + Ae is equal to the original matrix. | |
void | EliminateBC (const Array< int > &ess_dofs, DiagonalPolicy diag_policy) |
Eliminate essential (Dirichlet) boundary conditions. | |
void | SetDiagIdentity () |
If a row contains only one diag entry of zero, set it to 1. | |
virtual void | EliminateZeroRows (const real_t threshold=1e-12) |
If a row contains only zeros, set its diagonal to 1. | |
void | Gauss_Seidel_forw (const Vector &x, Vector &y) const |
Gauss-Seidel forward and backward iterations over a vector x. | |
void | Gauss_Seidel_back (const Vector &x, Vector &y) const |
real_t | GetJacobiScaling () const |
Determine appropriate scaling for Jacobi iteration. | |
void | Jacobi (const Vector &b, const Vector &x0, Vector &x1, real_t sc, bool use_abs_diag=false) const |
void | DiagScale (const Vector &b, Vector &x, real_t sc=1.0, bool use_abs_diag=false) const |
x = sc b / A_ii. When use_abs_diag = true, |A_ii| is used. | |
void | Jacobi2 (const Vector &b, const Vector &x0, Vector &x1, real_t sc=1.0) const |
void | Jacobi3 (const Vector &b, const Vector &x0, Vector &x1, real_t sc=1.0) const |
virtual void | Finalize (int skip_zeros=1) |
Finalize the matrix initialization, switching the storage format from LIL to CSR. | |
void | Finalize (int skip_zeros, bool fix_empty_rows) |
A slightly more general version of the Finalize(int) method. | |
bool | Finalized () const |
Returns whether or not CSR format has been finalized. | |
bool | ColumnsAreSorted () const |
Returns whether or not the columns are sorted. | |
void | Threshold (real_t tol, bool fix_empty_rows=false) |
Remove entries smaller in absolute value than a given tolerance tol. If fix_empty_rows is true, a zero value is inserted in the diagonal entry (for square matrices only) | |
void | GetBlocks (Array2D< SparseMatrix * > &blocks) const |
void | GetSubMatrix (const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const |
void | SetColPtr (const int row) const |
Initialize the SparseMatrix for fast access to the entries of the given row which becomes the "current row". | |
void | ClearColPtr () const |
Reset the "current row" set by calling SetColPtr(). This method must be called between any two calls to SetColPtr(). | |
real_t & | SearchRow (const int col) |
Perform a fast search for an entry in the "current row". See SetColPtr(). | |
void | _Add_ (const int col, const real_t a) |
Add a value to an entry in the "current row". See SetColPtr(). | |
void | _Set_ (const int col, const real_t a) |
Set an entry in the "current row". See SetColPtr(). | |
real_t | _Get_ (const int col) const |
Read the value of an entry in the "current row". See SetColPtr(). | |
real_t & | SearchRow (const int row, const int col) |
void | _Add_ (const int row, const int col, const real_t a) |
void | _Set_ (const int row, const int col, const real_t a) |
void | Set (const int i, const int j, const real_t val) |
void | Add (const int i, const int j, const real_t val) |
void | SetSubMatrix (const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1) |
void | SetSubMatrixTranspose (const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1) |
void | AddSubMatrix (const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1) |
bool | RowIsEmpty (const int row) const |
virtual int | GetRow (const int row, Array< int > &cols, Vector &srow) const |
Extract all column indices and values from a given row. | |
void | SetRow (const int row, const Array< int > &cols, const Vector &srow) |
void | AddRow (const int row, const Array< int > &cols, const Vector &srow) |
void | ScaleRow (const int row, const real_t scale) |
void | ScaleRows (const Vector &sl) |
this = diag(sl) * this; | |
void | ScaleColumns (const Vector &sr) |
this = this * diag(sr); | |
SparseMatrix & | operator+= (const SparseMatrix &B) |
Add the sparse matrix 'B' to '*this'. This operation will cause an error if '*this' is finalized and 'B' has larger sparsity pattern. | |
void | Add (const real_t a, const SparseMatrix &B) |
Add the sparse matrix 'B' scaled by the scalar 'a' into '*this'. Only entries in the sparsity pattern of '*this' are added. | |
SparseMatrix & | operator= (real_t a) |
SparseMatrix & | operator*= (real_t a) |
void | Print (std::ostream &out=mfem::out, int width_=4) const |
Prints matrix to stream out. | |
virtual void | PrintMatlab (std::ostream &out=mfem::out) const |
Prints matrix in matlab format. | |
void | PrintMM (std::ostream &out=mfem::out) const |
Prints matrix in Matrix Market sparse format. | |
void | PrintCSR (std::ostream &out) const |
Prints matrix to stream out in hypre_CSRMatrix format. | |
void | PrintCSR2 (std::ostream &out) const |
Prints a sparse matrix to stream out in CSR format. | |
void | PrintInfo (std::ostream &out) const |
Print various sparse matrix statistics. | |
real_t | IsSymmetric () const |
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix. | |
void | Symmetrize () |
(*this) = 1/2 ((*this) + (*this)^t) | |
virtual int | NumNonZeroElems () const |
Returns the number of the nonzero elements in the matrix. | |
real_t | MaxNorm () const |
int | CountSmallElems (real_t tol) const |
Count the number of entries with |a_ij| <= tol. | |
int | CheckFinite () const |
Count the number of entries that are NOT finite, i.e. Inf or Nan. | |
void | SetGraphOwner (bool ownij) |
Set the graph ownership flag (I and J arrays). | |
void | SetDataOwner (bool owna) |
Set the data ownership flag (A array). | |
bool | OwnsGraph () const |
Get the graph ownership flag (I and J arrays). | |
bool | OwnsData () const |
Get the data ownership flag (A array). | |
void | LoseData () |
Lose the ownership of the graph (I, J) and data (A) arrays. | |
void | Swap (SparseMatrix &other) |
virtual | ~SparseMatrix () |
Destroys sparse matrix. | |
Type | GetType () const |
Public Member Functions inherited from mfem::AbstractSparseMatrix | |
AbstractSparseMatrix (int s=0) | |
Creates a square matrix of the given size. | |
AbstractSparseMatrix (int h, int w) | |
Creates a matrix of the given height and width. | |
virtual | ~AbstractSparseMatrix () |
Destroys AbstractSparseMatrix. | |
Public Member Functions inherited from mfem::Matrix | |
Matrix (int s) | |
Creates a square matrix of size s. | |
Matrix (int h, int w) | |
Creates a matrix of the given height and width. | |
bool | IsSquare () const |
Returns whether the matrix is a square matrix. | |
virtual | ~Matrix () |
Destroys matrix. | |
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 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 | ~Operator () |
Virtual destructor. | |
Type | GetType () const |
Return the type ID of the Operator class. | |
Protected Types | |
typedef MemAlloc< RowNode, 1024 > | RowNodeAlloc |
Protected Member Functions | |
void | Destroy () |
void | SetEmpty () |
void | InitGPUSparse () |
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 | |
RowNode ** | Rows |
Array of linked lists, one for every row. This array represents the linked list (LIL) storage format. | |
int | current_row |
int * | ColPtrJ |
RowNode ** | ColPtrNode |
SparseMatrix * | At |
Transpose of A. Owned. Used to perform MultTranspose() on devices. | |
RowNodeAlloc * | NodesMem |
bool | isSorted |
Are the columns sorted already. | |
bool | useGPUSparse = true |
bool | initBuffers = false |
cusparseStatus_t | status |
cusparseMatDescr_t | descr = 0 |
cusparseSpMatDescr_t | matA_descr |
cusparseDnVecDescr_t | vecX_descr |
cusparseDnVecDescr_t | vecY_descr |
cusparseMatDescr_t | matA_descr |
hipsparseStatus_t | status |
hipsparseMatDescr_t | descr = 0 |
hipsparseSpMatDescr_t | matA_descr |
hipsparseDnVecDescr_t | vecX_descr |
hipsparseDnVecDescr_t | vecY_descr |
Arrays used by the CSR storage format. | |
Memory< int > | I |
Array with size (height+1) containing the row offsets. | |
Memory< int > | J |
Array with size I[height], containing the column indices for all matrix entries, as indexed by the I array. | |
Memory< real_t > | A |
Array with size I[height], containing the actual entries of the sparse matrix, as indexed by the I array. | |
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. | |
Static Protected Attributes | |
static int | SparseMatrixCount = 0 |
static size_t | bufferSize |
static void * | dBuffer = nullptr |
static cusparseHandle_t | handle |
static hipsparseHandle_t | handle |
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... | |
Data type sparse matrix.
Definition at line 50 of file sparsemat.hpp.
|
protected |
Definition at line 83 of file sparsemat.hpp.
|
inline |
Create an empty SparseMatrix.
Definition at line 131 of file sparsemat.hpp.
|
explicit |
Create a sparse matrix with flexible sparsity structure using a row-wise linked list (LIL) format.
New entries are added as needed by methods like AddSubMatrix(), SetSubMatrix(), etc. Calling Finalize() will convert the SparseMatrix to the more compact compressed sparse row (CSR) format.
Definition at line 98 of file sparsemat.cpp.
mfem::SparseMatrix::SparseMatrix | ( | int * | i, |
int * | j, | ||
real_t * | data, | ||
int | m, | ||
int | n ) |
Create a sparse matrix in CSR format. Ownership of i, j, and data is transferred to the SparseMatrix.
Definition at line 124 of file sparsemat.cpp.
mfem::SparseMatrix::SparseMatrix | ( | int * | i, |
int * | j, | ||
real_t * | data, | ||
int | m, | ||
int | n, | ||
bool | ownij, | ||
bool | owna, | ||
bool | issorted ) |
Create a sparse matrix in CSR format. Ownership of i, j, and data is optionally transferred to the SparseMatrix.
If the parameter data is NULL, then the internal A array is allocated by this constructor (initializing it with zeros and taking ownership, regardless of the parameter owna).
Definition at line 143 of file sparsemat.cpp.
mfem::SparseMatrix::SparseMatrix | ( | int | nrows, |
int | ncols, | ||
int | rowsize ) |
mfem::SparseMatrix::SparseMatrix | ( | const SparseMatrix & | mat, |
bool | copy_graph = true, | ||
MemoryType | mt = MemoryType::PRESERVE ) |
Copy constructor (deep copy).
If mat is finalized and copy_graph is false, the I and J arrays will use a shallow copy (copy the pointers only) without transferring ownership. If mt is MemoryType::PRESERVE the memory type of the resulting SparseMatrix's I, J, and A arrays will be the same as mat, otherwise the type will be mt for those arrays that are deep copied.
Definition at line 199 of file sparsemat.cpp.
mfem::SparseMatrix::SparseMatrix | ( | const Vector & | v | ) |
Create a SparseMatrix with diagonal v, i.e. A = Diag(v)
Definition at line 268 of file sparsemat.cpp.
|
virtual |
Destroys sparse matrix.
Definition at line 4303 of file sparsemat.cpp.
|
inline |
Add a value to an entry in the "current row". See SetColPtr().
Definition at line 585 of file sparsemat.hpp.
|
inline |
Definition at line 594 of file sparsemat.hpp.
|
inline |
Read the value of an entry in the "current row". See SetColPtr().
Definition at line 878 of file sparsemat.hpp.
|
inline |
Set an entry in the "current row". See SetColPtr().
Definition at line 588 of file sparsemat.hpp.
|
inline |
Definition at line 596 of file sparsemat.hpp.
y = |A| * x, using entry-wise absolute values of matrix A
Definition at line 1133 of file sparsemat.cpp.
y = |At| * x, using entry-wise absolute values of the transpose of matrix A
If the matrix is modified, call ResetTranspose() and optionally EnsureMultTranspose() to make sure this method uses the correct updated transpose.
Definition at line 1182 of file sparsemat.cpp.
int mfem::SparseMatrix::ActualWidth | ( | ) | const |
Returns the actual Width of the matrix.
This method can be called for matrices finalized or not.
Definition at line 3570 of file sparsemat.cpp.
void mfem::SparseMatrix::Add | ( | const int | i, |
const int | j, | ||
const real_t | val ) |
Definition at line 2842 of file sparsemat.cpp.
void mfem::SparseMatrix::Add | ( | const real_t | a, |
const SparseMatrix & | B ) |
Add the sparse matrix 'B' scaled by the scalar 'a' into '*this'. Only entries in the sparsity pattern of '*this' are added.
Definition at line 3226 of file sparsemat.cpp.
|
virtual |
y += A * x (default) or y += a * A * x
Implements mfem::AbstractSparseMatrix.
Definition at line 762 of file sparsemat.cpp.
|
virtual |
y += At * x (default) or y += a * At * x
If the matrix is modified, call ResetTranspose() and optionally EnsureMultTranspose() to make sure this method uses the correct updated transpose.
Implements mfem::AbstractSparseMatrix.
Definition at line 961 of file sparsemat.cpp.
Definition at line 3077 of file sparsemat.cpp.
void mfem::SparseMatrix::AddSubMatrix | ( | const Array< int > & | rows, |
const Array< int > & | cols, | ||
const DenseMatrix & | subm, | ||
int | skip_zeros = 1 ) |
Insert the DenseMatrix into this SparseMatrix at the specified rows and columns. If skip_zeros==0
, all entries from the DenseMatrix are added including zeros. If skip_zeros==2
, no zeros are added to the SparseMatrix regardless of their position in the matrix. Otherwise, the default skip_zeros
behavior is to omit the zero from the SparseMatrix unless it would break the symmetric structure of the SparseMatrix.
Definition at line 2777 of file sparsemat.cpp.
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
The actual values stored in the data array, A, are not used - this means that all entries in the sparsity pattern are considered to be true by this method.
Definition at line 1080 of file sparsemat.cpp.
y = At * x, treating all entries as booleans (zero=false, nonzero=true).
The actual values stored in the data array, A, are not used - this means that all entries in the sparsity pattern are considered to be true by this method.
Definition at line 1110 of file sparsemat.cpp.
void mfem::SparseMatrix::BuildTranspose | ( | ) | const |
Build and store internally the transpose of this matrix which will be used in the methods AddMultTranspose(), MultTranspose(), and AbsMultTranspose().
If this method has been called, the internal transpose matrix will be used to perform the action of the transpose matrix in AddMultTranspose(), MultTranspose(), and AbsMultTranspose().
Warning: any changes in this matrix will invalidate the internal transpose. To rebuild the transpose, call ResetTranspose() followed by (optionally) a call to this method. If the internal transpose is already built, this method has no effect.
When any non-serial-CPU backend is enabled, i.e. the call Device::Allows(~ Backend::CPU_MASK) returns true, the above methods require the internal transpose to be built. If that is not the case (i.e. the internal transpose is not built), these methods will automatically call EnsureMultTranspose(). When using any backend from Backend::CPU_MASK, calling this method is optional.
This method can only be used when the sparse matrix is finalized.
Definition at line 1014 of file sparsemat.cpp.
int mfem::SparseMatrix::CheckFinite | ( | ) | const |
Count the number of entries that are NOT finite, i.e. Inf or Nan.
Definition at line 1718 of file sparsemat.cpp.
|
inline |
Clear the contents of the SparseMatrix.
Definition at line 209 of file sparsemat.hpp.
|
inline |
Reset the "current row" set by calling SetColPtr(). This method must be called between any two calls to SetColPtr().
Definition at line 832 of file sparsemat.hpp.
|
inline |
Deprecated equivalent of ClearGPUSparse().
Definition at line 216 of file sparsemat.hpp.
void mfem::SparseMatrix::ClearGPUSparse | ( | ) |
Clear the cuSPARSE/hipSPARSE descriptors. This must be called after releasing the device memory of A.
Definition at line 81 of file sparsemat.cpp.
|
inline |
Returns whether or not the columns are sorted.
Definition at line 554 of file sparsemat.hpp.
int mfem::SparseMatrix::CountSmallElems | ( | real_t | tol | ) | const |
Count the number of entries with |a_ij| <= tol.
Definition at line 1690 of file sparsemat.cpp.
|
protected |
Definition at line 3537 of file sparsemat.cpp.
void mfem::SparseMatrix::DiagScale | ( | const Vector & | b, |
Vector & | x, | ||
real_t | sc = 1.0, | ||
bool | use_abs_diag = false ) const |
x = sc b / A_ii. When use_abs_diag = true, |A_ii| is used.
Definition at line 2675 of file sparsemat.cpp.
|
virtual |
|
virtual |
Returns constant reference to a_{ij}.
Implements mfem::Matrix.
Definition at line 628 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateBC | ( | const Array< int > & | ess_dofs, |
DiagonalPolicy | diag_policy ) |
Eliminate essential (Dirichlet) boundary conditions.
[in] | ess_dofs | indices of the degrees of freedom belonging to the essential boundary conditions. |
[in] | diag_policy | policy for diagonal entries. |
Definition at line 2370 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateCol | ( | int | col, |
DiagonalPolicy | dpolicy = DIAG_ZERO ) |
Eliminates the column col from the matrix.
Definition at line 1795 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateCols | ( | const Array< int > & | col_marker, |
SparseMatrix & | Ae ) |
Similar to EliminateCols + save the eliminated entries into Ae so that (*this) + Ae is equal to the original matrix.
Definition at line 1875 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateCols | ( | const Array< int > & | cols, |
const Vector * | x = NULL, | ||
Vector * | b = NULL ) |
Eliminate all columns i for which cols[i] != 0.
Elimination of a column means that all entries in the column are set to zero. In addition, if the pointers x and b are not NULL, the eliminated matrix entries are multiplied by the corresponding solution value in *x and subtracted from the r.h.s. vector, *b.
Definition at line 1836 of file sparsemat.cpp.
Eliminates a column from the transpose matrix.
Definition at line 1747 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRow | ( | int | row, |
DiagonalPolicy | dpolicy = DIAG_ZERO ) |
Eliminates a row from the matrix.
Definition at line 1763 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRowCol | ( | int | rc, |
const real_t | sol, | ||
Vector & | rhs, | ||
DiagonalPolicy | dpolicy = DIAG_ONE ) |
Eliminate row rc and column rc and modify the rhs using sol.
Eliminates the column rc to the rhs, deletes the row rc and replaces the element (rc,rc) with 1.0; assumes that element (i,rc) is assembled if and only if the element (rc,i) is assembled. By default, elements (rc,rc) are set to 1.0, although this behavior can be adjusted by changing the dpolicy parameter.
Definition at line 1909 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRowCol | ( | int | rc, |
DiagonalPolicy | dpolicy = DIAG_ONE ) |
Eliminate row rc and column rc.
Definition at line 2132 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRowCol | ( | int | rc, |
SparseMatrix & | Ae, | ||
DiagonalPolicy | dpolicy = DIAG_ONE ) |
Similar to EliminateRowCol(int, DiagonalPolicy) + save the eliminated entries into Ae so that (*this) + Ae is equal to the original matrix.
Definition at line 2276 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRowColDiag | ( | int | rc, |
real_t | value ) |
Perform elimination and set the diagonal entry to the given value.
Definition at line 2213 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRowColMultipleRHS | ( | int | rc, |
const Vector & | sol, | ||
DenseMatrix & | rhs, | ||
DiagonalPolicy | dpolicy = DIAG_ONE ) |
Similar to EliminateRowCol(int, const double, Vector &, DiagonalPolicy), but multiple values for eliminated unknowns are accepted, and accordingly multiple right-hand-sides are used.
Definition at line 2008 of file sparsemat.cpp.
|
virtual |
If a row contains only zeros, set its diagonal to 1.
Implements mfem::AbstractSparseMatrix.
Definition at line 2424 of file sparsemat.cpp.
|
inline |
Check if the SparseMatrix is empty.
Definition at line 219 of file sparsemat.hpp.
void mfem::SparseMatrix::EnsureMultTranspose | ( | ) | const |
Ensures that the matrix is capable of performing MultTranspose(), AddMultTranspose(), and AbsMultTranspose().
For non-serial-CPU backends (e.g. GPU, OpenMP), multiplying by the transpose requires that the internal transpose matrix be already built. When such a backend is enabled, this function will build the internal transpose matrix, see BuildTranspose().
For the serial CPU backends, the internal transpose is not required, and this function is a no-op. This allows for significant memory savings when the internal transpose matrix is not required.
Definition at line 1028 of file sparsemat.cpp.
void mfem::SparseMatrix::Finalize | ( | int | skip_zeros, |
bool | fix_empty_rows ) |
A slightly more general version of the Finalize(int) method.
Definition at line 1376 of file sparsemat.cpp.
|
inlinevirtual |
Finalize the matrix initialization, switching the storage format from LIL to CSR.
This method should be called once, after the matrix has been initialized. Internally, this method converts the matrix from row-wise linked list (LIL) format into CSR (compressed sparse row) format.
Reimplemented from mfem::Matrix.
Definition at line 546 of file sparsemat.hpp.
|
inline |
Returns whether or not CSR format has been finalized.
Definition at line 552 of file sparsemat.hpp.
Definition at line 2527 of file sparsemat.cpp.
Gauss-Seidel forward and backward iterations over a vector x.
Definition at line 2443 of file sparsemat.cpp.
void mfem::SparseMatrix::GetBlocks | ( | Array2D< SparseMatrix * > & | blocks | ) | const |
Split the matrix into M x N blocks of sparse matrices in CSR format. The 'blocks' array is M x N (i.e. M and N are determined by its dimensions) and its entries are overwritten by the new blocks.
Definition at line 1491 of file sparsemat.cpp.
|
inline |
Return the element data, i.e. the array A.
Definition at line 232 of file sparsemat.hpp.
|
inline |
Return the element data, i.e. the array A, const version.
Definition at line 234 of file sparsemat.hpp.
void mfem::SparseMatrix::GetDiag | ( | Vector & | d | ) | const |
Returns the Diagonal of A.
Definition at line 691 of file sparsemat.cpp.
|
inline |
Return the array I.
Definition at line 222 of file sparsemat.hpp.
|
inline |
Return the array I, const version.
Definition at line 224 of file sparsemat.hpp.
|
inline |
Return the array J.
Definition at line 227 of file sparsemat.hpp.
|
inline |
Return the array J, const version.
Definition at line 229 of file sparsemat.hpp.
real_t mfem::SparseMatrix::GetJacobiScaling | ( | ) | const |
Determine appropriate scaling for Jacobi iteration.
Definition at line 2610 of file sparsemat.cpp.
|
inlinevirtual |
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.
Definition at line 342 of file sparsemat.hpp.
Definition at line 269 of file sparsemat.hpp.
Definition at line 270 of file sparsemat.hpp.
|
inline |
Definition at line 237 of file sparsemat.hpp.
|
inline |
Definition at line 238 of file sparsemat.hpp.
|
inline |
Definition at line 253 of file sparsemat.hpp.
|
inline |
Definition at line 254 of file sparsemat.hpp.
Extract all column indices and values from a given row.
If the matrix is finalized (i.e. in CSR format), cols and srow will simply be references to the specific portion of the J and A arrays. As required by the AbstractSparseMatrix interface this method returns:
Implements mfem::AbstractSparseMatrix.
Definition at line 2990 of file sparsemat.cpp.
int * mfem::SparseMatrix::GetRowColumns | ( | const int | row | ) |
Return a pointer to the column indices in a row.
Definition at line 391 of file sparsemat.cpp.
const int * mfem::SparseMatrix::GetRowColumns | ( | const int | row | ) | const |
Return a pointer to the column indices in a row, const version.
Definition at line 398 of file sparsemat.cpp.
real_t * mfem::SparseMatrix::GetRowEntries | ( | const int | row | ) |
Return a pointer to the entries in a row.
Definition at line 405 of file sparsemat.cpp.
const real_t * mfem::SparseMatrix::GetRowEntries | ( | const int | row | ) | const |
Return a pointer to the entries in a row, const version.
Definition at line 412 of file sparsemat.cpp.
real_t mfem::SparseMatrix::GetRowNorml1 | ( | int | irow | ) | const |
For i = irow compute \( x_i = \sum_j | A_{i, j} | \).
Definition at line 1291 of file sparsemat.cpp.
void mfem::SparseMatrix::GetRowSums | ( | Vector & | x | ) | const |
For all i compute \( x_i = \sum_j A_{ij} \).
Definition at line 1268 of file sparsemat.cpp.
void mfem::SparseMatrix::GetSubMatrix | ( | const Array< int > & | rows, |
const Array< int > & | cols, | ||
DenseMatrix & | subm ) const |
Definition at line 2941 of file sparsemat.cpp.
|
inline |
Definition at line 713 of file sparsemat.hpp.
|
inline |
Definition at line 277 of file sparsemat.hpp.
|
inline |
Definition at line 245 of file sparsemat.hpp.
|
inline |
Definition at line 261 of file sparsemat.hpp.
|
inline |
Definition at line 281 of file sparsemat.hpp.
|
inline |
Definition at line 249 of file sparsemat.hpp.
|
inline |
Definition at line 265 of file sparsemat.hpp.
|
inline |
Definition at line 279 of file sparsemat.hpp.
|
inline |
Definition at line 247 of file sparsemat.hpp.
|
inline |
Definition at line 263 of file sparsemat.hpp.
|
protected |
Definition at line 64 of file sparsemat.cpp.
Compute y^t A x.
Definition at line 1227 of file sparsemat.cpp.
|
virtual |
This virtual method is not supported: it always returns NULL.
Implements mfem::Matrix.
Definition at line 1742 of file sparsemat.cpp.
real_t mfem::SparseMatrix::IsSymmetric | ( | ) | const |
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix.
Definition at line 1581 of file sparsemat.cpp.
void mfem::SparseMatrix::Jacobi | ( | const Vector & | b, |
const Vector & | x0, | ||
Vector & | x1, | ||
real_t | sc, | ||
bool | use_abs_diag = false ) const |
One scaled Jacobi iteration for the system A x = b. x1 = x0 + sc D^{-1} (b - A x0) where D is the diag of A. Absolute values of D are used when use_abs_diag = true.
Definition at line 2643 of file sparsemat.cpp.
void mfem::SparseMatrix::Jacobi2 | ( | const Vector & | b, |
const Vector & | x0, | ||
Vector & | x1, | ||
real_t | sc = 1.0 ) const |
x1 = x0 + sc D^{-1} (b - A x0) where \( D_{ii} = \sum_j |A_{ij}| \).
Definition at line 2763 of file sparsemat.cpp.
void mfem::SparseMatrix::Jacobi3 | ( | const Vector & | b, |
const Vector & | x0, | ||
Vector & | x1, | ||
real_t | sc = 1.0 ) const |
x1 = x0 + sc D^{-1} (b - A x0) where \( D_{ii} = \sum_j A_{ij} \).
Definition at line 2770 of file sparsemat.cpp.
|
inline |
Lose the ownership of the graph (I, J) and data (A) arrays.
Definition at line 706 of file sparsemat.hpp.
void mfem::SparseMatrix::MakeRef | ( | const SparseMatrix & | master | ) |
Clear the contents of the SparseMatrix and make it a reference to master.
After this call, the matrix will point to the same data as master but it will not own its data. The master must be finalized.
Definition at line 313 of file sparsemat.cpp.
real_t mfem::SparseMatrix::MaxNorm | ( | ) | const |
Definition at line 1665 of file sparsemat.cpp.
int mfem::SparseMatrix::MaxRowSize | ( | ) | const |
Returns the maximum number of elements among all rows.
Definition at line 367 of file sparsemat.cpp.
void mfem::SparseMatrix::MoveDiagonalFirst | ( | ) |
Move the diagonal entry to the first position in each row, preserving the order of the rest of the columns.
Definition at line 599 of file sparsemat.cpp.
Matrix vector multiplication.
Implements mfem::AbstractSparseMatrix.
Definition at line 755 of file sparsemat.cpp.
Multiply a vector with the transposed matrix. y = At * x.
If the matrix is modified, call ResetTranspose() and optionally EnsureMultTranspose() to make sure this method uses the correct updated transpose.
Implements mfem::AbstractSparseMatrix.
Definition at line 954 of file sparsemat.cpp.
|
virtual |
Returns the number of the nonzero elements in the matrix.
Implements mfem::AbstractSparseMatrix.
Definition at line 1642 of file sparsemat.cpp.
real_t & mfem::SparseMatrix::operator() | ( | int | i, |
int | j ) |
Returns reference to A[i][j].
Definition at line 633 of file sparsemat.cpp.
const real_t & mfem::SparseMatrix::operator() | ( | int | i, |
int | j ) const |
Returns reference to A[i][j].
Definition at line 656 of file sparsemat.cpp.
SparseMatrix & mfem::SparseMatrix::operator*= | ( | real_t | a | ) |
Definition at line 3275 of file sparsemat.cpp.
SparseMatrix & mfem::SparseMatrix::operator+= | ( | const SparseMatrix & | B | ) |
Add the sparse matrix 'B' to '*this'. This operation will cause an error if '*this' is finalized and 'B' has larger sparsity pattern.
Definition at line 3196 of file sparsemat.cpp.
SparseMatrix & mfem::SparseMatrix::operator= | ( | const SparseMatrix & | rhs | ) |
Assignment operator: deep copy.
Definition at line 303 of file sparsemat.cpp.
SparseMatrix & mfem::SparseMatrix::operator= | ( | real_t | a | ) |
Definition at line 3249 of file sparsemat.cpp.
void mfem::SparseMatrix::OverrideSize | ( | int | height_, |
int | width_ ) |
Sets the height and width of the matrix.
This function should generally be called when manually constructing the CSR I, J, and A arrays in conjunction with the SparseMatrix::SparseMatrix() constructor.
Definition at line 297 of file sparsemat.cpp.
|
inline |
Get the data ownership flag (A array).
Definition at line 703 of file sparsemat.hpp.
|
inline |
Get the graph ownership flag (I and J arrays).
Definition at line 700 of file sparsemat.hpp.
void mfem::SparseMatrix::PartAddMult | ( | const Array< int > & | rows, |
const Vector & | x, | ||
Vector & | y, | ||
const real_t | a = 1.0 ) const |
Definition at line 1062 of file sparsemat.cpp.
Definition at line 1036 of file sparsemat.cpp.
|
virtual |
Prints matrix to stream out.
Reimplemented from mfem::Matrix.
Definition at line 3299 of file sparsemat.cpp.
void mfem::SparseMatrix::PrintCSR | ( | std::ostream & | out | ) | const |
Prints matrix to stream out in hypre_CSRMatrix format.
Definition at line 3429 of file sparsemat.cpp.
void mfem::SparseMatrix::PrintCSR2 | ( | std::ostream & | out | ) | const |
Prints a sparse matrix to stream out in CSR format.
Definition at line 3457 of file sparsemat.cpp.
void mfem::SparseMatrix::PrintInfo | ( | std::ostream & | out | ) | const |
Print various sparse matrix statistics.
Definition at line 3486 of file sparsemat.cpp.
|
virtual |
Prints matrix in matlab format.
Reimplemented from mfem::Operator.
Definition at line 3347 of file sparsemat.cpp.
void mfem::SparseMatrix::PrintMM | ( | std::ostream & | out = mfem::out | ) | const |
Prints matrix in Matrix Market sparse format.
Definition at line 3388 of file sparsemat.cpp.
|
inline |
Definition at line 271 of file sparsemat.hpp.
|
inline |
Definition at line 239 of file sparsemat.hpp.
|
inline |
Definition at line 255 of file sparsemat.hpp.
|
inline |
Definition at line 275 of file sparsemat.hpp.
|
inline |
Definition at line 243 of file sparsemat.hpp.
|
inline |
Definition at line 259 of file sparsemat.hpp.
void mfem::SparseMatrix::ResetTranspose | ( | ) | const |
Reset (destroy) the internal transpose matrix. See BuildTranspose() for more details.
Definition at line 1022 of file sparsemat.cpp.
bool mfem::SparseMatrix::RowIsEmpty | ( | const int | row | ) | const |
Definition at line 2969 of file sparsemat.cpp.
int mfem::SparseMatrix::RowSize | ( | const int | i | ) | const |
Returns the number of elements in row i.
Definition at line 344 of file sparsemat.cpp.
void mfem::SparseMatrix::ScaleColumns | ( | const Vector & | sr | ) |
this = this * diag(sr);
Definition at line 3168 of file sparsemat.cpp.
void mfem::SparseMatrix::ScaleRow | ( | const int | row, |
const real_t | scale ) |
Definition at line 3109 of file sparsemat.cpp.
void mfem::SparseMatrix::ScaleRows | ( | const Vector & | sl | ) |
this = diag(sl) * this;
Definition at line 3137 of file sparsemat.cpp.
|
inline |
Perform a fast search for an entry in the "current row". See SetColPtr().
If the matrix is not finalized and the entry is not found in the SparseMatrix, it will be added to the sparsity pattern initialized with zero. If the matrix is finalized and the entry is not found, an error will be generated.
Definition at line 851 of file sparsemat.hpp.
|
inline |
Definition at line 892 of file sparsemat.hpp.
void mfem::SparseMatrix::Set | ( | const int | i, |
const int | j, | ||
const real_t | val ) |
Definition at line 2823 of file sparsemat.cpp.
|
inline |
Initialize the SparseMatrix for fast access to the entries of the given row which becomes the "current row".
Fast access to the entries of the "current row" can be performed using the methods: SearchRow(const int), Add(const int, const double), Set(const int, const double), and Get(const int).
Definition at line 797 of file sparsemat.hpp.
|
inline |
Set the data ownership flag (A array).
Definition at line 697 of file sparsemat.hpp.
void mfem::SparseMatrix::SetDiagIdentity | ( | ) |
If a row contains only one diag entry of zero, set it to 1.
Definition at line 2413 of file sparsemat.cpp.
|
protected |
Definition at line 325 of file sparsemat.cpp.
|
inline |
Set the graph ownership flag (I and J arrays).
Definition at line 693 of file sparsemat.hpp.
Definition at line 3030 of file sparsemat.cpp.
void mfem::SparseMatrix::SetSubMatrix | ( | const Array< int > & | rows, |
const Array< int > & | cols, | ||
const DenseMatrix & | subm, | ||
int | skip_zeros = 1 ) |
Definition at line 2861 of file sparsemat.cpp.
void mfem::SparseMatrix::SetSubMatrixTranspose | ( | const Array< int > & | rows, |
const Array< int > & | cols, | ||
const DenseMatrix & | subm, | ||
int | skip_zeros = 1 ) |
Definition at line 2900 of file sparsemat.cpp.
void mfem::SparseMatrix::SetWidth | ( | int | width_ = -1 | ) |
Change the width of a SparseMatrix.
If width_ = -1 (DEFAULT), this routine will set the new width to the actual Width of the matrix awidth = max(J) + 1. Values 0 <= width_ < awidth are not allowed (error check in Debug Mode only)
This method can be called for matrices finalized or not.
Definition at line 419 of file sparsemat.cpp.
|
inline |
For backward compatibility, define Size() to be synonym of Height().
Definition at line 206 of file sparsemat.hpp.
void mfem::SparseMatrix::SortColumnIndices | ( | ) |
Sort the column indices corresponding to each row.
Definition at line 457 of file sparsemat.cpp.
void mfem::SparseMatrix::Swap | ( | SparseMatrix & | other | ) |
Definition at line 4283 of file sparsemat.cpp.
void mfem::SparseMatrix::Symmetrize | ( | ) |
(*this) = 1/2 ((*this) + (*this)^t)
Definition at line 1623 of file sparsemat.cpp.
void mfem::SparseMatrix::Threshold | ( | real_t | tol, |
bool | fix_empty_rows = false ) |
Remove entries smaller in absolute value than a given tolerance tol. If fix_empty_rows is true, a zero value is inserted in the diagonal entry (for square matrices only)
Definition at line 1315 of file sparsemat.cpp.
DenseMatrix * mfem::SparseMatrix::ToDenseMatrix | ( | ) | const |
Produces a DenseMatrix from a SparseMatrix.
Definition at line 725 of file sparsemat.cpp.
void mfem::SparseMatrix::ToDenseMatrix | ( | DenseMatrix & | B | ) | const |
Produces a DenseMatrix from a SparseMatrix.
Definition at line 738 of file sparsemat.cpp.
|
inline |
Deprecated equivalent of UseGPUSparse().
Definition at line 194 of file sparsemat.hpp.
|
inline |
Runtime option to use cuSPARSE or hipSPARSE. Only valid when using a CUDA or HIP backend.
Definition at line 191 of file sparsemat.hpp.
|
inline |
Definition at line 273 of file sparsemat.hpp.
|
inline |
Definition at line 241 of file sparsemat.hpp.
|
inline |
Definition at line 257 of file sparsemat.hpp.
Array with size I[height], containing the actual entries of the sparse matrix, as indexed by the I array.
Definition at line 68 of file sparsemat.hpp.
|
mutableprotected |
Transpose of A. Owned. Used to perform MultTranspose() on devices.
Definition at line 80 of file sparsemat.hpp.
|
staticprotected |
Definition at line 101 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 76 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 77 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 75 of file sparsemat.hpp.
|
staticprotected |
Definition at line 102 of file sparsemat.hpp.
|
protected |
Definition at line 108 of file sparsemat.hpp.
|
protected |
Definition at line 121 of file sparsemat.hpp.
|
staticprotected |
Definition at line 107 of file sparsemat.hpp.
|
staticprotected |
Definition at line 120 of file sparsemat.hpp.
|
protected |
Array with size (height+1) containing the row offsets.
The data for row r, 0 <= r < height, is at offsets j, I[r] <= j < I[r+1]. The offsets, j, are indices in the J and A arrays. The first entry in this array is always zero, I[0] = 0, and the last entry, I[height], gives the total number of entries stored (at a minimum, all nonzeros must be represented) in the sparse matrix.
Definition at line 62 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 103 of file sparsemat.hpp.
|
protected |
Are the columns sorted already.
Definition at line 88 of file sparsemat.hpp.
|
protected |
Array with size I[height], containing the column indices for all matrix entries, as indexed by the I array.
Definition at line 65 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 111 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 115 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 123 of file sparsemat.hpp.
|
protected |
Definition at line 84 of file sparsemat.hpp.
|
protected |
Array of linked lists, one for every row. This array represents the linked list (LIL) storage format.
Definition at line 73 of file sparsemat.hpp.
|
staticprotected |
Definition at line 100 of file sparsemat.hpp.
|
protected |
Definition at line 106 of file sparsemat.hpp.
|
protected |
Definition at line 119 of file sparsemat.hpp.
|
protected |
Definition at line 93 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 112 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 124 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 113 of file sparsemat.hpp.
|
mutableprotected |
Definition at line 125 of file sparsemat.hpp.