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

Data type sparse matrix. More...

#include <sparsemat.hpp>

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

Public Member Functions

 SparseMatrix ()
 Create an empty SparseMatrix. More...
 
 SparseMatrix (int nrows, int ncols=-1)
 Create a sparse matrix with flexible sparsity structure using a row-wise linked list (LIL) format. More...
 
 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. More...
 
 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. More...
 
 SparseMatrix (int nrows, int ncols, int rowsize)
 Create a sparse matrix in CSR format where each row has space allocated for exactly rowsize entries. More...
 
 SparseMatrix (const SparseMatrix &mat, bool copy_graph=true)
 Copy constructor (deep copy). More...
 
SparseMatrixoperator= (const SparseMatrix &rhs)
 Assignment operator: deep copy. More...
 
void MakeRef (const SparseMatrix &master)
 Clear the contents of the SparseMatrix and make it a reference to master. More...
 
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 () const
 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, i.e. array A. 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...
 
void MoveDiagonalFirst ()
 Move the diagonal entry to the first position in each row, preserving the order of the rest of the columns. 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 \( x_i = \sum_j A_{ij} \). More...
 
double GetRowNorml1 (int irow) const
 For i = irow compute \( x_i = \sum_j | A_{i, j} | \). More...
 
virtual MatrixInverseInverse () const
 This virtual method is not supported: it always returns NULL. More...
 
void EliminateRow (int row, const double sol, Vector &rhs)
 Eliminates a column from the transpose matrix. More...
 
void EliminateRow (int row, DiagonalPolicy dpolicy=DIAG_ZERO)
 Eliminates a row from the matrix. More...
 
void EliminateCol (int col, DiagonalPolicy dpolicy=DIAG_ZERO)
 Eliminates the column col from the matrix. More...
 
void EliminateCols (const Array< int > &cols, const 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, DiagonalPolicy dpolicy=DIAG_ONE)
 Eliminate row rc and column rc and modify the rhs using sol. More...
 
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. More...
 
void EliminateRowColDiag (int rc, double value)
 Perform elimination and set the diagonal entry to the given value. More...
 
void EliminateRowCol (int rc, DiagonalPolicy dpolicy=DIAG_ONE)
 Eliminate row rc and column rc. More...
 
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. More...
 
void SetDiagIdentity ()
 If a row contains only one diag entry of zero, set it to 1. More...
 
virtual void EliminateZeroRows (const double threshold=1e-12)
 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)
 Finalize the matrix initialization, switching the storage format from LIL to CSR. More...
 
void Finalize (int skip_zeros, bool fix_empty_rows)
 A slightly more general version of the Finalize(int) method. More...
 
bool Finalized () const
 
bool areColumnsSorted () const
 
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
 
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
 Extract all column indices and values from a given row. More...
 
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)
 
SparseMatrixoperator+= (const SparseMatrix &B)
 
void Add (const double a, const SparseMatrix &B)
 
SparseMatrixoperator= (double a)
 
SparseMatrixoperator*= (double a)
 
void Print (std::ostream &out=mfem::out, int width_=4) const
 Prints matrix to stream out. More...
 
void PrintMatlab (std::ostream &out=mfem::out) const
 Prints matrix in matlab format. More...
 
void PrintMM (std::ostream &out=mfem::out) 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...
 
void PrintInfo (std::ostream &out) const
 Print various sparse matrix staticstics. 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...
 
int CheckFinite () const
 Count the number of entries that are NOT finite, i.e. Inf or Nan. 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...
 
Type GetType () const
 
- 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)
 Construct an Operator with the given height (output size) and width (input size). More...
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows(). More...
 
int NumRows () const
 Get the number of rows (size of output) of the Operator. Synonym with Height(). More...
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols(). More...
 
int NumCols () const
 Get the number of columns (size of input) of the Operator. Synonym with Width(). More...
 
virtual OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error. More...
 
virtual const OperatorGetProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More...
 
virtual const OperatorGetRestriction () const
 Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More...
 
void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
 Form a constrained linear system using a matrix-free approach. More...
 
virtual void RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x)
 Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear system obtained from Operator::FormLinearSystem(). More...
 
void PrintMatlab (std::ostream &out, int n=0, int m=0) const
 Prints operator with input size n and output size m in Matlab format. More...
 
virtual ~Operator ()
 Virtual destructor. More...
 
Type GetType () const
 Return the type ID of the Operator class. More...
 

Protected Types

typedef MemAlloc< RowNode, 1024 > RowNodeAlloc
 

Protected Member Functions

void Destroy ()
 
void SetEmpty ()
 

Protected Attributes

RowNode ** Rows
 Array of linked lists, one for every row. This array represents the linked list (LIL) storage format. More...
 
int current_row
 
int * ColPtrJ
 
RowNode ** ColPtrNode
 
RowNodeAllocNodesMem
 
bool ownGraph
 Say whether we own the pointers for I and J (should we free them?). More...
 
bool ownData
 Say whether we own the pointers for A (should we free them?). More...
 
bool isSorted
 Are the columns sorted already. More...
 
Arrays used by the CSR storage format.
int * I
 Array with size (height+1) containing the row offsets. More...
 
int * J
 Array with size I[height], containing the column indices for all matrix entries, as indexed by the I array. More...
 
double * A
 Array with size I[height], containing the actual entries of the sparse matrix, as indexed by the I array. More...
 
- Protected Attributes inherited from mfem::Operator
int height
 Dimension of the output / number of rows in the matrix. More...
 
int width
 Dimension of the input / number of columns in the matrix. More...
 

Additional Inherited Members

- Public Types inherited from mfem::Matrix
enum  DiagonalPolicy { DIAG_ZERO, DIAG_ONE, DIAG_KEEP }
 
- Public Types inherited from mfem::Operator
enum  Type {
  ANY_TYPE, MFEM_SPARSEMAT, Hypre_ParCSR, PETSC_MATAIJ,
  PETSC_MATIS, PETSC_MATSHELL, PETSC_MATNEST, PETSC_MATHYPRE,
  PETSC_MATGENERIC
}
 Enumeration defining IDs for some classes derived from Operator. More...
 

Detailed Description

Data type sparse matrix.

Definition at line 38 of file sparsemat.hpp.

Member Typedef Documentation

typedef MemAlloc<RowNode, 1024> mfem::SparseMatrix::RowNodeAlloc
protected

Definition at line 68 of file sparsemat.hpp.

Constructor & Destructor Documentation

mfem::SparseMatrix::SparseMatrix ( )
inline

Create an empty SparseMatrix.

Definition at line 85 of file sparsemat.hpp.

mfem::SparseMatrix::SparseMatrix ( int  nrows,
int  ncols = -1 
)
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 30 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 53 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 70 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 99 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 121 of file sparsemat.cpp.

virtual mfem::SparseMatrix::~SparseMatrix ( )
inlinevirtual

Destroys sparse matrix.

Definition at line 433 of file sparsemat.hpp.

Member Function Documentation

void mfem::SparseMatrix::_Add_ ( const int  col,
const double  a 
)
inline

Definition at line 324 of file sparsemat.hpp.

void mfem::SparseMatrix::_Add_ ( const int  row,
const int  col,
const double  a 
)
inline

Definition at line 331 of file sparsemat.hpp.

double mfem::SparseMatrix::_Get_ ( const int  col) const
inline

Definition at line 569 of file sparsemat.hpp.

void mfem::SparseMatrix::_Set_ ( const int  col,
const double  a 
)
inline

Definition at line 326 of file sparsemat.hpp.

void mfem::SparseMatrix::_Set_ ( const int  row,
const int  col,
const double  a 
)
inline

Definition at line 333 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 2696 of file sparsemat.cpp.

void mfem::SparseMatrix::Add ( const int  i,
const int  j,
const double  a 
)

Definition at line 2035 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 2407 of file sparsemat.cpp.

void mfem::SparseMatrix::AddMult ( const Vector x,
Vector y,
const double  a = 1.0 
) const
virtual

y += A * x (default) or y += a * A * x

Implements mfem::AbstractSparseMatrix.

Definition at line 492 of file sparsemat.cpp.

void mfem::SparseMatrix::AddMultTranspose ( const Vector x,
Vector y,
const double  a = 1.0 
) const
virtual

y += At * x (default) or y += a * At * x

Implements mfem::AbstractSparseMatrix.

Definition at line 569 of file sparsemat.cpp.

void mfem::SparseMatrix::AddRow ( const int  row,
const Array< int > &  cols,
const Vector srow 
)

Definition at line 2258 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 1978 of file sparsemat.cpp.

bool mfem::SparseMatrix::areColumnsSorted ( ) const
inline

Definition at line 311 of file sparsemat.hpp.

void mfem::SparseMatrix::BooleanMult ( const Array< int > &  x,
Array< int > &  y 
) const

y = A * x, but treat all elements as booleans (zero=false, nonzero=true).

Definition at line 644 of file sparsemat.cpp.

void mfem::SparseMatrix::BooleanMultTranspose ( const Array< int > &  x,
Array< int > &  y 
) const

y = At * x, but treat all elements as booleans (zero=false, nonzero=true).

Definition at line 667 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 1056 of file sparsemat.cpp.

void mfem::SparseMatrix::Clear ( )
inline

Clear the contents of the SparseMatrix.

Definition at line 128 of file sparsemat.hpp.

void mfem::SparseMatrix::ClearColPtr ( ) const
inline

Definition at line 527 of file sparsemat.hpp.

int mfem::SparseMatrix::CountSmallElems ( double  tol) const

Count the number of entries with |a_ij| <= tol.

Definition at line 1028 of file sparsemat.cpp.

void mfem::SparseMatrix::Destroy ( )
protected

Definition at line 2648 of file sparsemat.cpp.

void mfem::SparseMatrix::DiagScale ( const Vector b,
Vector x,
double  sc = 1.0 
) const

Definition at line 1898 of file sparsemat.cpp.

double & mfem::SparseMatrix::Elem ( int  i,
int  j 
)
virtual

Returns reference to a_{ij}.

Implements mfem::Matrix.

Definition at line 390 of file sparsemat.cpp.

const double & mfem::SparseMatrix::Elem ( int  i,
int  j 
) const
virtual

Returns constant reference to a_{ij}.

Implements mfem::Matrix.

Definition at line 395 of file sparsemat.cpp.

void mfem::SparseMatrix::EliminateCol ( int  col,
DiagonalPolicy  dpolicy = DIAG_ZERO 
)

Eliminates the column col from the matrix.

  • If dpolicy = DIAG_ZERO, all entries in the column will be set to 0.
  • If dpolicy = DIAG_ONE (matrix must be square), the diagonal entry will be set equal to 1 and all other entries in the column to 0.
  • The policy DIAG_KEEP is not supported.

Definition at line 1133 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 1173 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 1085 of file sparsemat.cpp.

void mfem::SparseMatrix::EliminateRow ( int  row,
DiagonalPolicy  dpolicy = DIAG_ZERO 
)

Eliminates a row from the matrix.

  • If dpolicy = DIAG_ZERO, all the entries in the row will be set to 0.
  • If dpolicy = DIAG_ONE (matrix must be square), the diagonal entry will be set equal to 1 and all other entries in the row to 0.
  • The policy DIAG_KEEP is not supported.

Definition at line 1101 of file sparsemat.cpp.

void mfem::SparseMatrix::EliminateRowCol ( int  rc,
const double  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 1204 of file sparsemat.cpp.

void mfem::SparseMatrix::EliminateRowCol ( int  rc,
DiagonalPolicy  dpolicy = DIAG_ONE 
)

Eliminate row rc and column rc.

Definition at line 1424 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 1554 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 1497 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 1300 of file sparsemat.cpp.

void mfem::SparseMatrix::EliminateZeroRows ( const double  threshold = 1e-12)
virtual

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

Implements mfem::AbstractSparseMatrix.

Definition at line 1657 of file sparsemat.cpp.

bool mfem::SparseMatrix::Empty ( ) const
inline

Check if the SparseMatrix is empty.

Definition at line 131 of file sparsemat.hpp.

virtual void mfem::SparseMatrix::Finalize ( int  skip_zeros = 1)
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 305 of file sparsemat.hpp.

void mfem::SparseMatrix::Finalize ( int  skip_zeros,
bool  fix_empty_rows 
)

A slightly more general version of the Finalize(int) method.

Definition at line 755 of file sparsemat.cpp.

bool mfem::SparseMatrix::Finalized ( ) const
inline

Definition at line 310 of file sparsemat.hpp.

void mfem::SparseMatrix::Gauss_Seidel_back ( const Vector x,
Vector y 
) const

Definition at line 1759 of file sparsemat.cpp.

void mfem::SparseMatrix::Gauss_Seidel_forw ( const Vector x,
Vector y 
) const

Gauss-Seidel forward and backward iterations over a vector x.

Definition at line 1684 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 836 of file sparsemat.cpp.

double* mfem::SparseMatrix::GetData ( ) const
inline

Return element data, i.e. array A.

Definition at line 138 of file sparsemat.hpp.

void mfem::SparseMatrix::GetDiag ( Vector d) const

Returns the Diagonal of A.

Definition at line 458 of file sparsemat.cpp.

int* mfem::SparseMatrix::GetI ( ) const
inline

Return the array I.

Definition at line 134 of file sparsemat.hpp.

int* mfem::SparseMatrix::GetJ ( ) const
inline

Return the array J.

Definition at line 136 of file sparsemat.hpp.

double mfem::SparseMatrix::GetJacobiScaling ( ) const

Determine appropriate scaling for Jacobi iteration.

Definition at line 1834 of file sparsemat.cpp.

int mfem::SparseMatrix::GetRow ( const int  row,
Array< int > &  cols,
Vector srow 
) const
virtual

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 2171 of file sparsemat.cpp.

int * mfem::SparseMatrix::GetRowColumns ( const int  row)

Return a pointer to the column indices in a row.

Definition at line 271 of file sparsemat.cpp.

const int * mfem::SparseMatrix::GetRowColumns ( const int  row) const

Definition at line 278 of file sparsemat.cpp.

double * mfem::SparseMatrix::GetRowEntries ( const int  row)

Return a pointer to the entries in a row.

Definition at line 285 of file sparsemat.cpp.

const double * mfem::SparseMatrix::GetRowEntries ( const int  row) const

Definition at line 292 of file sparsemat.cpp.

double mfem::SparseMatrix::GetRowNorml1 ( int  irow) const

For i = irow compute \( x_i = \sum_j | A_{i, j} | \).

Definition at line 735 of file sparsemat.cpp.

void mfem::SparseMatrix::GetRowSums ( Vector x) const

For all i compute \( x_i = \sum_j A_{ij} \).

Definition at line 716 of file sparsemat.cpp.

void mfem::SparseMatrix::GetSubMatrix ( const Array< int > &  rows,
const Array< int > &  cols,
DenseMatrix subm 
) const

Definition at line 2122 of file sparsemat.cpp.

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

Definition at line 435 of file sparsemat.hpp.

double mfem::SparseMatrix::InnerProduct ( const Vector x,
const Vector y 
) const

Compute y^t A x.

Definition at line 690 of file sparsemat.cpp.

MatrixInverse * mfem::SparseMatrix::Inverse ( ) const
virtual

This virtual method is not supported: it always returns NULL.

Implements mfem::Matrix.

Definition at line 1080 of file sparsemat.cpp.

double mfem::SparseMatrix::IsSymmetric ( ) const

Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix.

Definition at line 926 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 1867 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 \( D_{ii} = \sum_j |A_{ij}| \).

Definition at line 1930 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 \( D_{ii} = \sum_j A_{ij} \).

Definition at line 1954 of file sparsemat.cpp.

void mfem::SparseMatrix::LoseData ( )
inline

Lose the ownership of the graph (I, J) and data (A) arrays.

Definition at line 428 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 197 of file sparsemat.cpp.

double mfem::SparseMatrix::MaxNorm ( ) const

Definition at line 1005 of file sparsemat.cpp.

int mfem::SparseMatrix::MaxRowSize ( ) const

Returns the maximum number of elements among all rows.

Definition at line 247 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 366 of file sparsemat.cpp.

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

Matrix vector multiplication.

Implements mfem::AbstractSparseMatrix.

Definition at line 486 of file sparsemat.cpp.

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

Multiply a vector with the transposed matrix. y = At * x.

Implements mfem::AbstractSparseMatrix.

Definition at line 563 of file sparsemat.cpp.

int mfem::SparseMatrix::NumNonZeroElems ( ) const
virtual

Returns the number of the nonzero elements in the matrix.

Implements mfem::AbstractSparseMatrix.

Definition at line 983 of file sparsemat.cpp.

double & mfem::SparseMatrix::operator() ( int  i,
int  j 
)

Returns reference to A[i][j].

Definition at line 400 of file sparsemat.cpp.

const double & mfem::SparseMatrix::operator() ( int  i,
int  j 
) const

Returns reference to A[i][j].

Definition at line 423 of file sparsemat.cpp.

SparseMatrix & mfem::SparseMatrix::operator*= ( double  a)

Definition at line 2448 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 2377 of file sparsemat.cpp.

SparseMatrix & mfem::SparseMatrix::operator= ( const SparseMatrix rhs)

Assignment operator: deep copy.

Definition at line 187 of file sparsemat.cpp.

SparseMatrix & mfem::SparseMatrix::operator= ( double  a)

Definition at line 2430 of file sparsemat.cpp.

bool mfem::SparseMatrix::OwnsData ( ) const
inline

Get the data ownership flag (A array).

Definition at line 426 of file sparsemat.hpp.

bool mfem::SparseMatrix::OwnsGraph ( ) const
inline

Get the graph ownership flag (I and J arrays).

Definition at line 424 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 626 of file sparsemat.cpp.

void mfem::SparseMatrix::PartMult ( const Array< int > &  rows,
const Vector x,
Vector y 
) const

Definition at line 608 of file sparsemat.cpp.

void mfem::SparseMatrix::Print ( std::ostream &  out = mfem::out,
int  width_ = 4 
) const
virtual

Prints matrix to stream out.

Reimplemented from mfem::Matrix.

Definition at line 2466 of file sparsemat.cpp.

void mfem::SparseMatrix::PrintCSR ( std::ostream &  out) const

Prints matrix to stream out in hypre_CSRMatrix format.

Definition at line 2548 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 2572 of file sparsemat.cpp.

void mfem::SparseMatrix::PrintInfo ( std::ostream &  out) const

Print various sparse matrix staticstics.

Definition at line 2597 of file sparsemat.cpp.

void mfem::SparseMatrix::PrintMatlab ( std::ostream &  out = mfem::out) const

Prints matrix in matlab format.

Definition at line 2510 of file sparsemat.cpp.

void mfem::SparseMatrix::PrintMM ( std::ostream &  out = mfem::out) const

Prints matrix in Matrix Market sparse format.

Definition at line 2528 of file sparsemat.cpp.

bool mfem::SparseMatrix::RowIsEmpty ( const int  row) const

Definition at line 2150 of file sparsemat.cpp.

int mfem::SparseMatrix::RowSize ( const int  i) const

Returns the number of elements in row i.

Definition at line 224 of file sparsemat.cpp.

void mfem::SparseMatrix::ScaleColumns ( const Vector sr)

Definition at line 2349 of file sparsemat.cpp.

void mfem::SparseMatrix::ScaleRow ( const int  row,
const double  scale 
)

Definition at line 2290 of file sparsemat.cpp.

void mfem::SparseMatrix::ScaleRows ( const Vector sl)

Definition at line 2318 of file sparsemat.cpp.

double & mfem::SparseMatrix::SearchRow ( const int  col)
inline

Definition at line 542 of file sparsemat.hpp.

double & mfem::SparseMatrix::SearchRow ( const int  row,
const int  col 
)
inline

Definition at line 583 of file sparsemat.hpp.

void mfem::SparseMatrix::Set ( const int  i,
const int  j,
const double  a 
)

Definition at line 2016 of file sparsemat.cpp.

void mfem::SparseMatrix::SetColPtr ( const int  row) const
inline

Definition at line 492 of file sparsemat.hpp.

void mfem::SparseMatrix::SetDataOwner ( bool  owna)
inline

Set the data ownership flag (A array).

Definition at line 422 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 1648 of file sparsemat.cpp.

void mfem::SparseMatrix::SetEmpty ( )
protected

Definition at line 209 of file sparsemat.cpp.

void mfem::SparseMatrix::SetGraphOwner ( bool  ownij)
inline

Set the graph ownership flag (I and J arrays).

Definition at line 420 of file sparsemat.hpp.

void mfem::SparseMatrix::SetRow ( const int  row,
const Array< int > &  cols,
const Vector srow 
)

Definition at line 2210 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 2054 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 2087 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 299 of file sparsemat.cpp.

int mfem::SparseMatrix::Size ( ) const
inline

For backward compatibility, define Size() to be synonym of Height().

Definition at line 125 of file sparsemat.hpp.

void mfem::SparseMatrix::SortColumnIndices ( )

Sort the column indices corresponding to each row.

Definition at line 337 of file sparsemat.cpp.

void mfem::SparseMatrix::Swap ( SparseMatrix other)

Definition at line 3260 of file sparsemat.cpp.

void mfem::SparseMatrix::Symmetrize ( )

(*this) = 1/2 ((*this) + (*this)^t)

Definition at line 968 of file sparsemat.cpp.

int mfem::SparseMatrix::Walk ( int &  i,
int &  j,
double &  a 
)

Walks the sparse matrix.

Member Data Documentation

double* mfem::SparseMatrix::A
protected

Array with size I[height], containing the actual entries of the sparse matrix, as indexed by the I array.

Definition at line 56 of file sparsemat.hpp.

int* mfem::SparseMatrix::ColPtrJ
mutableprotected

Definition at line 64 of file sparsemat.hpp.

RowNode** mfem::SparseMatrix::ColPtrNode
mutableprotected

Definition at line 65 of file sparsemat.hpp.

int mfem::SparseMatrix::current_row
mutableprotected

Definition at line 63 of file sparsemat.hpp.

int* mfem::SparseMatrix::I
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 50 of file sparsemat.hpp.

bool mfem::SparseMatrix::isSorted
protected

Are the columns sorted already.

Definition at line 78 of file sparsemat.hpp.

int* mfem::SparseMatrix::J
protected

Array with size I[height], containing the column indices for all matrix entries, as indexed by the I array.

Definition at line 53 of file sparsemat.hpp.

RowNodeAlloc* mfem::SparseMatrix::NodesMem
protected

Definition at line 69 of file sparsemat.hpp.

bool mfem::SparseMatrix::ownData
protected

Say whether we own the pointers for A (should we free them?).

Definition at line 75 of file sparsemat.hpp.

bool mfem::SparseMatrix::ownGraph
protected

Say whether we own the pointers for I and J (should we free them?).

Definition at line 73 of file sparsemat.hpp.

RowNode** mfem::SparseMatrix::Rows
protected

Array of linked lists, one for every row. This array represents the linked list (LIL) storage format.

Definition at line 61 of file sparsemat.hpp.


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