MFEM
v3.2
Finite element discretization library
|
Data type sparse matrix. More...
#include <sparsemat.hpp>
Public Member Functions | |
SparseMatrix () | |
Create an empty SparseMatrix. More... | |
SparseMatrix (int nrows, int ncols=-1) | |
SparseMatrix (int *i, int *j, double *data, int m, int n) | |
SparseMatrix (int *i, int *j, double *data, int m, int n, bool ownij, bool owna, bool issorted) | |
SparseMatrix (int nrows, int ncols, int rowsize) | |
SparseMatrix (const SparseMatrix &mat, bool copy_graph=true) | |
void | MakeRef (const SparseMatrix &master) |
int | Size () const |
For backward compatibility define Size to be synonym of Height() More... | |
void | Clear () |
Clear the contents of the SparseMatrix. More... | |
bool | Empty () |
Check if the SparseMatrix is empty. More... | |
int * | GetI () const |
Return the array I. More... | |
int * | GetJ () const |
Return the array J. More... | |
double * | GetData () const |
Return element data. More... | |
int | RowSize (const int i) const |
Returns the number of elements in row i. More... | |
int | MaxRowSize () const |
Returns the maximum number of elements among all rows. More... | |
int * | GetRowColumns (const int row) |
Return a pointer to the column indices in a row. More... | |
const int * | GetRowColumns (const int row) const |
double * | GetRowEntries (const int row) |
Return a pointer to the entries in a row. More... | |
const double * | GetRowEntries (const int row) const |
void | SetWidth (int width_=-1) |
Change the width of a SparseMatrix. More... | |
int | ActualWidth () |
Returns the actual Width of the matrix. More... | |
void | SortColumnIndices () |
Sort the column indices corresponding to each row. More... | |
virtual double & | Elem (int i, int j) |
Returns reference to a_{ij}. More... | |
virtual const double & | Elem (int i, int j) const |
Returns constant reference to a_{ij}. More... | |
double & | operator() (int i, int j) |
Returns reference to A[i][j]. More... | |
const double & | operator() (int i, int j) const |
Returns reference to A[i][j]. More... | |
void | GetDiag (Vector &d) const |
Returns the Diagonal of A. More... | |
virtual void | Mult (const Vector &x, Vector &y) const |
Matrix vector multiplication. More... | |
void | AddMult (const Vector &x, Vector &y, const double a=1.0) const |
y += A * x (default) or y += a * A * x More... | |
void | MultTranspose (const Vector &x, Vector &y) const |
Multiply a vector with the transposed matrix. y = At * x. More... | |
void | AddMultTranspose (const Vector &x, Vector &y, const double a=1.0) const |
y += At * x (default) or y += a * At * x More... | |
void | PartMult (const Array< int > &rows, const Vector &x, Vector &y) const |
void | PartAddMult (const Array< int > &rows, const Vector &x, Vector &y, const double a=1.0) const |
void | BooleanMult (const Array< int > &x, Array< int > &y) const |
y = A * x, but treat all elements as booleans (zero=false, nonzero=true). More... | |
void | BooleanMultTranspose (const Array< int > &x, Array< int > &y) const |
y = At * x, but treat all elements as booleans (zero=false, nonzero=true). More... | |
double | InnerProduct (const Vector &x, const Vector &y) const |
Compute y^t A x. More... | |
void | GetRowSums (Vector &x) const |
For all i compute . More... | |
double | GetRowNorml1 (int irow) const |
For i = irow compute . More... | |
virtual MatrixInverse * | Inverse () const |
Returns a pointer to approximation of the matrix inverse. More... | |
void | EliminateRow (int row, const double sol, Vector &rhs) |
Eliminates a column from the transpose matrix. More... | |
void | EliminateRow (int row, int setOneDiagonal=0) |
Eliminates a row from the matrix. More... | |
void | EliminateCol (int col) |
void | EliminateCols (Array< int > &cols, Vector *x=NULL, Vector *b=NULL) |
Eliminate all columns 'i' for which cols[i] != 0. More... | |
void | EliminateRowCol (int rc, const double sol, Vector &rhs, int d=0) |
void | EliminateRowColMultipleRHS (int rc, const Vector &sol, DenseMatrix &rhs, int d=0) |
void | EliminateRowCol (int rc, int d=0) |
void | EliminateRowColDiag (int rc, double value) |
Perform elimination and set the diagonal entry to the given value. More... | |
void | EliminateRowCol (int rc, SparseMatrix &Ae, int d=0) |
void | SetDiagIdentity () |
If a row contains only one diag entry of zero, set it to 1. More... | |
void | EliminateZeroRows () |
If a row contains only zeros, set its diagonal to 1. More... | |
void | Gauss_Seidel_forw (const Vector &x, Vector &y) const |
Gauss-Seidel forward and backward iterations over a vector x. More... | |
void | Gauss_Seidel_back (const Vector &x, Vector &y) const |
double | GetJacobiScaling () const |
Determine appropriate scaling for Jacobi iteration. More... | |
void | Jacobi (const Vector &b, const Vector &x0, Vector &x1, double sc) const |
void | DiagScale (const Vector &b, Vector &x, double sc=1.0) const |
void | Jacobi2 (const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const |
void | Jacobi3 (const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const |
virtual void | Finalize (int skip_zeros=1) |
bool | Finalized () const |
bool | areColumnsSorted () const |
void | GetBlocks (Array2D< SparseMatrix * > &blocks) const |
void | GetSubMatrix (const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) |
void | SetColPtr (const int row) const |
void | ClearColPtr () const |
double & | SearchRow (const int col) |
void | _Add_ (const int col, const double a) |
void | _Set_ (const int col, const double a) |
double | _Get_ (const int col) const |
double & | SearchRow (const int row, const int col) |
void | _Add_ (const int row, const int col, const double a) |
void | _Set_ (const int row, const int col, const double a) |
void | Set (const int i, const int j, const double a) |
void | Add (const int i, const int j, const double a) |
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 |
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 double scale) |
void | ScaleRows (const Vector &sl) |
void | ScaleColumns (const Vector &sr) |
SparseMatrix & | operator+= (const SparseMatrix &B) |
void | Add (const double a, const SparseMatrix &B) |
SparseMatrix & | operator= (double a) |
SparseMatrix & | operator*= (double a) |
void | Print (std::ostream &out=std::cout, int width_=4) const |
Prints matrix to stream out. More... | |
void | PrintMatlab (std::ostream &out=std::cout) const |
Prints matrix in matlab format. More... | |
void | PrintMM (std::ostream &out=std::cout) const |
Prints matrix in Matrix Market sparse format. More... | |
void | PrintCSR (std::ostream &out) const |
Prints matrix to stream out in hypre_CSRMatrix format. More... | |
void | PrintCSR2 (std::ostream &out) const |
Prints a sparse matrix to stream out in CSR format. More... | |
int | Walk (int &i, int &j, double &a) |
Walks the sparse matrix. More... | |
double | IsSymmetric () const |
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix. More... | |
void | Symmetrize () |
(*this) = 1/2 ((*this) + (*this)^t) More... | |
virtual int | NumNonZeroElems () const |
Returns the number of the nonzero elements in the matrix. More... | |
double | MaxNorm () const |
int | CountSmallElems (double tol) const |
Count the number of entries with |a_ij| < tol. More... | |
void | SetGraphOwner (bool ownij) |
Set the graph ownership flag (I and J arrays). More... | |
void | SetDataOwner (bool owna) |
Set the data ownership flag (A array). More... | |
bool | OwnsGraph () const |
Get the graph ownership flag (I and J arrays). More... | |
bool | OwnsData () const |
Get the data ownership flag (A array). More... | |
void | LoseData () |
Lose the ownership of the graph (I, J) and data (A) arrays. More... | |
void | Swap (SparseMatrix &other) |
virtual | ~SparseMatrix () |
Destroys sparse matrix. More... | |
Public Member Functions inherited from mfem::AbstractSparseMatrix | |
AbstractSparseMatrix (int s=0) | |
Creates a square matrix of the given size. More... | |
AbstractSparseMatrix (int h, int w) | |
Creates a matrix of the given height and width. More... | |
virtual | ~AbstractSparseMatrix () |
Destroys AbstractSparseMatrix. More... | |
Public Member Functions inherited from mfem::Matrix | |
Matrix (int s) | |
Creates a square matrix of size s. More... | |
Matrix (int h, int w) | |
Creates a matrix of the given height and width. More... | |
virtual | ~Matrix () |
Destroys matrix. More... | |
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) | |
int | Height () const |
Get the height (size of output) of the Operator. Synonym with NumRows. More... | |
int | NumRows () const |
int | Width () const |
Get the width (size of input) of the Operator. Synonym with NumCols. More... | |
int | NumCols () const |
virtual Operator & | GetGradient (const Vector &x) const |
Evaluate the gradient operator at the point x. 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 () |
Additional Inherited Members | |
Protected Attributes inherited from mfem::Operator | |
int | height |
int | width |
Data type sparse matrix.
Definition at line 38 of file sparsemat.hpp.
|
inline |
Create an empty SparseMatrix.
Definition at line 75 of file sparsemat.hpp.
|
explicit |
Create a sparse matrix with flexible sparsity structure using a row-wise linked list 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 29 of file sparsemat.cpp.
mfem::SparseMatrix::SparseMatrix | ( | int * | i, |
int * | j, | ||
double * | 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 52 of file sparsemat.cpp.
mfem::SparseMatrix::SparseMatrix | ( | int * | i, |
int * | j, | ||
double * | 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.
Definition at line 69 of file sparsemat.cpp.
mfem::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. SetRow can then be called or the I, J, A arrays can be used directly.
Definition at line 98 of file sparsemat.cpp.
mfem::SparseMatrix::SparseMatrix | ( | const SparseMatrix & | mat, |
bool | copy_graph = true |
||
) |
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.
Definition at line 120 of file sparsemat.cpp.
|
inlinevirtual |
Destroys sparse matrix.
Definition at line 378 of file sparsemat.hpp.
|
inline |
Definition at line 275 of file sparsemat.hpp.
|
inline |
Definition at line 282 of file sparsemat.hpp.
|
inline |
Definition at line 512 of file sparsemat.hpp.
|
inline |
Definition at line 277 of file sparsemat.hpp.
|
inline |
Definition at line 284 of file sparsemat.hpp.
int mfem::SparseMatrix::ActualWidth | ( | ) |
Returns the actual Width of the matrix.
This method can be called for matrices finalized or not.
Definition at line 2426 of file sparsemat.cpp.
void mfem::SparseMatrix::Add | ( | const int | i, |
const int | j, | ||
const double | a | ||
) |
Definition at line 1816 of file sparsemat.cpp.
void mfem::SparseMatrix::Add | ( | const double | 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 2188 of file sparsemat.cpp.
|
virtual |
y += A * x (default) or y += a * A * x
Implements mfem::AbstractSparseMatrix.
Definition at line 445 of file sparsemat.cpp.
|
virtual |
y += At * x (default) or y += a * At * x
Implements mfem::AbstractSparseMatrix.
Definition at line 520 of file sparsemat.cpp.
Definition at line 2039 of file sparsemat.cpp.
void mfem::SparseMatrix::AddSubMatrix | ( | const Array< int > & | rows, |
const Array< int > & | cols, | ||
const DenseMatrix & | subm, | ||
int | skip_zeros = 1 |
||
) |
Definition at line 1759 of file sparsemat.cpp.
|
inline |
Definition at line 262 of file sparsemat.hpp.
y = A * x, but treat all elements as booleans (zero=false, nonzero=true).
Definition at line 595 of file sparsemat.cpp.
y = At * x, but treat all elements as booleans (zero=false, nonzero=true).
Definition at line 618 of file sparsemat.cpp.
|
inline |
Clear the contents of the SparseMatrix.
Definition at line 111 of file sparsemat.hpp.
|
inline |
Definition at line 470 of file sparsemat.hpp.
int mfem::SparseMatrix::CountSmallElems | ( | double | tol | ) | const |
Count the number of entries with |a_ij| < tol.
Definition at line 948 of file sparsemat.cpp.
Definition at line 1679 of file sparsemat.cpp.
|
virtual |
|
virtual |
Returns constant reference to a_{ij}.
Implements mfem::Matrix.
Definition at line 359 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateCol | ( | int | col | ) |
Definition at line 1026 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateCols | ( | Array< int > & | cols, |
Vector * | x = NULL , |
||
Vector * | b = NULL |
||
) |
Eliminate all columns 'i' for which cols[i] != 0.
Definition at line 1040 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRow | ( | int | row, |
const double | sol, | ||
Vector & | rhs | ||
) |
Eliminates a column from the transpose matrix.
Definition at line 983 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRow | ( | int | row, |
int | setOneDiagonal = 0 |
||
) |
Eliminates a row from the matrix.
If setOneDiagonal = 0, all the entries in the row will be set to 0. If setOneDiagonal = 1 (matrix must be square), the diagonal entry will be set equal to 1 and all the others entries to 0.
Definition at line 999 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRowCol | ( | int | rc, |
const double | sol, | ||
Vector & | rhs, | ||
int | d = 0 |
||
) |
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. If d != 0 then the element (rc,rc) remains the same.
Definition at line 1071 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRowCol | ( | int | rc, |
int | d = 0 |
||
) |
Definition at line 1236 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRowCol | ( | int | rc, |
SparseMatrix & | Ae, | ||
int | d = 0 |
||
) |
Definition at line 1358 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRowColDiag | ( | int | rc, |
double | value | ||
) |
Perform elimination and set the diagonal entry to the given value.
Definition at line 1301 of file sparsemat.cpp.
void mfem::SparseMatrix::EliminateRowColMultipleRHS | ( | int | rc, |
const Vector & | sol, | ||
DenseMatrix & | rhs, | ||
int | d = 0 |
||
) |
Like previous one, but multiple values for eliminated dofs are accepted, and accordingly multiple right-hand-sides are used.
Definition at line 1151 of file sparsemat.cpp.
|
virtual |
If a row contains only zeros, set its diagonal to 1.
Implements mfem::AbstractSparseMatrix.
Definition at line 1438 of file sparsemat.cpp.
|
inline |
Check if the SparseMatrix is empty.
Definition at line 114 of file sparsemat.hpp.
|
virtual |
Finalize the matrix initialization. The function should be called only once, after the matrix has been initialized. Internally, this method converts the matrix from row-wise linked list format into CSR (compressed sparse row) format.
Reimplemented from mfem::Matrix.
Definition at line 706 of file sparsemat.cpp.
|
inline |
Definition at line 261 of file sparsemat.hpp.
Definition at line 1540 of file sparsemat.cpp.
Gauss-Seidel forward and backward iterations over a vector x.
Definition at line 1465 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 778 of file sparsemat.cpp.
|
inline |
Return element data.
Definition at line 121 of file sparsemat.hpp.
void mfem::SparseMatrix::GetDiag | ( | Vector & | d | ) | const |
Returns the Diagonal of A.
Definition at line 411 of file sparsemat.cpp.
|
inline |
Return the array I.
Definition at line 117 of file sparsemat.hpp.
|
inline |
Return the array J.
Definition at line 119 of file sparsemat.hpp.
double mfem::SparseMatrix::GetJacobiScaling | ( | ) | const |
Determine appropriate scaling for Jacobi iteration.
Definition at line 1615 of file sparsemat.cpp.
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: 0, if cols and srow are copies of the values in the matrix, i.e. when the matrix is open. 1, if cols and srow are views of the values in the matrix, i.e. when the matrix is finalized.
Implements mfem::AbstractSparseMatrix.
Definition at line 1952 of file sparsemat.cpp.
int * mfem::SparseMatrix::GetRowColumns | ( | const int | row | ) |
Return a pointer to the column indices in a row.
Definition at line 259 of file sparsemat.cpp.
const int * mfem::SparseMatrix::GetRowColumns | ( | const int | row | ) | const |
Definition at line 266 of file sparsemat.cpp.
double * mfem::SparseMatrix::GetRowEntries | ( | const int | row | ) |
Return a pointer to the entries in a row.
Definition at line 273 of file sparsemat.cpp.
const double * mfem::SparseMatrix::GetRowEntries | ( | const int | row | ) | const |
Definition at line 280 of file sparsemat.cpp.
double mfem::SparseMatrix::GetRowNorml1 | ( | int | irow | ) | const |
For i = irow compute .
Definition at line 686 of file sparsemat.cpp.
void mfem::SparseMatrix::GetRowSums | ( | Vector & | x | ) | const |
For all i compute .
Definition at line 667 of file sparsemat.cpp.
void mfem::SparseMatrix::GetSubMatrix | ( | const Array< int > & | rows, |
const Array< int > & | cols, | ||
DenseMatrix & | subm | ||
) |
Definition at line 1903 of file sparsemat.cpp.
Compute y^t A x.
Definition at line 641 of file sparsemat.cpp.
|
virtual |
Returns a pointer to approximation of the matrix inverse.
Implements mfem::Matrix.
Definition at line 978 of file sparsemat.cpp.
double mfem::SparseMatrix::IsSymmetric | ( | ) | const |
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix.
Definition at line 868 of file sparsemat.cpp.
void mfem::SparseMatrix::Jacobi | ( | const Vector & | b, |
const Vector & | x0, | ||
Vector & | x1, | ||
double | sc | ||
) | 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.
Definition at line 1648 of file sparsemat.cpp.
void mfem::SparseMatrix::Jacobi2 | ( | const Vector & | b, |
const Vector & | x0, | ||
Vector & | x1, | ||
double | sc = 1.0 |
||
) | const |
x1 = x0 + sc D^{-1} (b - A x0) where .
Definition at line 1711 of file sparsemat.cpp.
void mfem::SparseMatrix::Jacobi3 | ( | const Vector & | b, |
const Vector & | x0, | ||
Vector & | x1, | ||
double | sc = 1.0 |
||
) | const |
x1 = x0 + sc D^{-1} (b - A x0) where .
Definition at line 1735 of file sparsemat.cpp.
|
inline |
Lose the ownership of the graph (I, J) and data (A) arrays.
Definition at line 373 of file sparsemat.hpp.
void mfem::SparseMatrix::MakeRef | ( | const SparseMatrix & | master | ) |
Clear the contents of the SparseMatrix and make it a reference to 'master', i.e. it will point to the same data as 'master' but it will not own its data. The 'master' must be finalized.
Definition at line 186 of file sparsemat.cpp.
double mfem::SparseMatrix::MaxNorm | ( | ) | const |
Definition at line 925 of file sparsemat.cpp.
int mfem::SparseMatrix::MaxRowSize | ( | ) | const |
Returns the maximum number of elements among all rows.
Definition at line 235 of file sparsemat.cpp.
Matrix vector multiplication.
Implements mfem::AbstractSparseMatrix.
Definition at line 439 of file sparsemat.cpp.
Multiply a vector with the transposed matrix. y = At * x.
Implements mfem::AbstractSparseMatrix.
Definition at line 514 of file sparsemat.cpp.
|
virtual |
Returns the number of the nonzero elements in the matrix.
Implements mfem::AbstractSparseMatrix.
Definition at line 905 of file sparsemat.cpp.
double & mfem::SparseMatrix::operator() | ( | int | i, |
int | j | ||
) |
Returns reference to A[i][j].
Definition at line 364 of file sparsemat.cpp.
const double & mfem::SparseMatrix::operator() | ( | int | i, |
int | j | ||
) | const |
Returns reference to A[i][j].
Definition at line 388 of file sparsemat.cpp.
SparseMatrix & mfem::SparseMatrix::operator*= | ( | double | a | ) |
Definition at line 2229 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 2158 of file sparsemat.cpp.
SparseMatrix & mfem::SparseMatrix::operator= | ( | double | a | ) |
Definition at line 2211 of file sparsemat.cpp.
|
inline |
Get the data ownership flag (A array).
Definition at line 371 of file sparsemat.hpp.
|
inline |
Get the graph ownership flag (I and J arrays).
Definition at line 369 of file sparsemat.hpp.
void mfem::SparseMatrix::PartAddMult | ( | const Array< int > & | rows, |
const Vector & | x, | ||
Vector & | y, | ||
const double | a = 1.0 |
||
) | const |
Definition at line 577 of file sparsemat.cpp.
Definition at line 559 of file sparsemat.cpp.
|
virtual |
Prints matrix to stream out.
Reimplemented from mfem::Matrix.
Definition at line 2247 of file sparsemat.cpp.
void mfem::SparseMatrix::PrintCSR | ( | std::ostream & | out | ) | const |
Prints matrix to stream out in hypre_CSRMatrix format.
Definition at line 2329 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 2353 of file sparsemat.cpp.
void mfem::SparseMatrix::PrintMatlab | ( | std::ostream & | out = std::cout | ) | const |
Prints matrix in matlab format.
Definition at line 2291 of file sparsemat.cpp.
void mfem::SparseMatrix::PrintMM | ( | std::ostream & | out = std::cout | ) | const |
Prints matrix in Matrix Market sparse format.
Definition at line 2309 of file sparsemat.cpp.
bool mfem::SparseMatrix::RowIsEmpty | ( | const int | row | ) | const |
Definition at line 1931 of file sparsemat.cpp.
int mfem::SparseMatrix::RowSize | ( | const int | i | ) | const |
Returns the number of elements in row i.
Definition at line 212 of file sparsemat.cpp.
void mfem::SparseMatrix::ScaleColumns | ( | const Vector & | sr | ) |
Definition at line 2130 of file sparsemat.cpp.
void mfem::SparseMatrix::ScaleRow | ( | const int | row, |
const double | scale | ||
) |
Definition at line 2071 of file sparsemat.cpp.
void mfem::SparseMatrix::ScaleRows | ( | const Vector & | sl | ) |
Definition at line 2099 of file sparsemat.cpp.
|
inline |
Definition at line 485 of file sparsemat.hpp.
|
inline |
Definition at line 526 of file sparsemat.hpp.
void mfem::SparseMatrix::Set | ( | const int | i, |
const int | j, | ||
const double | a | ||
) |
Definition at line 1797 of file sparsemat.cpp.
|
inline |
Definition at line 435 of file sparsemat.hpp.
|
inline |
Set the data ownership flag (A array).
Definition at line 367 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 1429 of file sparsemat.cpp.
|
inline |
Set the graph ownership flag (I and J arrays).
Definition at line 365 of file sparsemat.hpp.
Definition at line 1991 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 1835 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 1868 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 287 of file sparsemat.cpp.
|
inline |
For backward compatibility define Size to be synonym of Height()
Definition at line 108 of file sparsemat.hpp.
void mfem::SparseMatrix::SortColumnIndices | ( | ) |
Sort the column indices corresponding to each row.
Definition at line 325 of file sparsemat.cpp.
void mfem::SparseMatrix::Swap | ( | SparseMatrix & | other | ) |
Definition at line 2990 of file sparsemat.cpp.
void mfem::SparseMatrix::Symmetrize | ( | ) |
(*this) = 1/2 ((*this) + (*this)^t)
Definition at line 890 of file sparsemat.cpp.
int mfem::SparseMatrix::Walk | ( | int & | i, |
int & | j, | ||
double & | a | ||
) |
Walks the sparse matrix.