MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Friends | 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 (int nrows, int ncols=0)
 Creates sparse matrix. More...
 
 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)
 
int Size () const
 For backward compatibility define Size to be synonym of Height() 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
 
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
 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 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 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)
 
SparseMatrixoperator+= (SparseMatrix &B)
 
void Add (const double a, const SparseMatrix &B)
 
SparseMatrixoperator= (double a)
 
SparseMatrixoperator*= (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 LoseData ()
 Call this if data has been stolen. More...
 
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 OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. More...
 
void PrintMatlab (std::ostream &out, int n=0, int m=0)
 Prints operator with input size n and output size m in matlab format. More...
 
virtual ~Operator ()
 

Friends

void Swap (SparseMatrix &A, SparseMatrix &B)
 

Additional Inherited Members

- Protected Attributes inherited from mfem::Operator
int height
 
int width
 

Detailed Description

Data type sparse matrix.

Definition at line 38 of file sparsemat.hpp.

Constructor & Destructor Documentation

mfem::SparseMatrix::SparseMatrix ( int  nrows,
int  ncols = 0 
)
explicit

Creates sparse matrix.

Definition at line 27 of file sparsemat.cpp.

mfem::SparseMatrix::SparseMatrix ( int *  i,
int *  j,
double *  data,
int  m,
int  n 
)

Definition at line 50 of file sparsemat.cpp.

mfem::SparseMatrix::SparseMatrix ( int *  i,
int *  j,
double *  data,
int  m,
int  n,
bool  ownij,
bool  owna,
bool  issorted 
)

Definition at line 67 of file sparsemat.cpp.

mfem::SparseMatrix::~SparseMatrix ( )
virtual

Destroys sparse matrix.

Definition at line 2067 of file sparsemat.cpp.

Member Function Documentation

int mfem::SparseMatrix::ActualWidth ( )

Returns the actual Width of the matrix.

This method can be called for matrices finalized or not.

Definition at line 2115 of file sparsemat.cpp.

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

Definition at line 1525 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 1877 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 318 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 393 of file sparsemat.cpp.

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

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

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

Definition at line 235 of file sparsemat.hpp.

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

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

Definition at line 736 of file sparsemat.cpp.

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

Definition at line 1388 of file sparsemat.cpp.

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

Returns reference to a_{ij}.

Implements mfem::Matrix.

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

void mfem::SparseMatrix::EliminateCol ( int  col)

Definition at line 814 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 828 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 771 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 787 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 859 of file sparsemat.cpp.

void mfem::SparseMatrix::EliminateRowCol ( int  rc,
int  d = 0 
)

Definition at line 1008 of file sparsemat.cpp.

void mfem::SparseMatrix::EliminateRowCol ( int  rc,
SparseMatrix Ae,
int  d = 0 
)

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

void mfem::SparseMatrix::EliminateZeroRows ( )
virtual

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

Implements mfem::AbstractSparseMatrix.

Definition at line 1147 of file sparsemat.cpp.

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

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

Definition at line 234 of file sparsemat.hpp.

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

Definition at line 1249 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 1174 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 601 of file sparsemat.cpp.

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

Return element data.

Definition at line 101 of file sparsemat.hpp.

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

Returns the Diagonal of A.

Definition at line 284 of file sparsemat.cpp.

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

Return the array I.

Definition at line 97 of file sparsemat.hpp.

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

Return the array J.

Definition at line 99 of file sparsemat.hpp.

double mfem::SparseMatrix::GetJacobiScaling ( ) const

Determine appropriate scaling for Jacobi iteration.

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

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

Return a pointer to the column indices in a row.

Definition at line 132 of file sparsemat.cpp.

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

Definition at line 139 of file sparsemat.cpp.

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

Return a pointer to the entries in a row.

Definition at line 146 of file sparsemat.cpp.

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

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

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

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

Definition at line 490 of file sparsemat.cpp.

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

Definition at line 1612 of file sparsemat.cpp.

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

Compute y^t A x.

Definition at line 468 of file sparsemat.cpp.

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

Returns a pointer to approximation of the matrix inverse.

Implements mfem::Matrix.

Definition at line 766 of file sparsemat.cpp.

double mfem::SparseMatrix::IsSymmetric ( ) const

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

Definition at line 656 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 1357 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 1420 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 1444 of file sparsemat.cpp.

void mfem::SparseMatrix::LoseData ( )
inline

Call this if data has been stolen.

Definition at line 322 of file sparsemat.hpp.

double mfem::SparseMatrix::MaxNorm ( ) const

Definition at line 713 of file sparsemat.cpp.

int mfem::SparseMatrix::MaxRowSize ( ) const

Returns the maximum number of elements among all rows.

Definition at line 108 of file sparsemat.cpp.

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

Matrix vector multiplication.

Implements mfem::AbstractSparseMatrix.

Definition at line 312 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 387 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 693 of file sparsemat.cpp.

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

Returns reference to A[i][j].

Definition at line 237 of file sparsemat.cpp.

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

Returns reference to A[i][j].

Definition at line 261 of file sparsemat.cpp.

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

Definition at line 1918 of file sparsemat.cpp.

SparseMatrix & mfem::SparseMatrix::operator+= ( 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 1847 of file sparsemat.cpp.

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

Definition at line 1900 of file sparsemat.cpp.

void mfem::SparseMatrix::PartAddMult ( const Array< int > &  rows,
const Vector x,
Vector y,
const double  a = 1.0 
) const

Definition at line 450 of file sparsemat.cpp.

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

Definition at line 432 of file sparsemat.cpp.

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

Prints matrix to stream out.

Reimplemented from mfem::Matrix.

Definition at line 1936 of file sparsemat.cpp.

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

Prints matrix to stream out in hypre_CSRMatrix format.

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

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

Prints matrix in matlab format.

Definition at line 1980 of file sparsemat.cpp.

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

Prints matrix in Matrix Market sparse format.

Definition at line 1998 of file sparsemat.cpp.

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

Definition at line 1640 of file sparsemat.cpp.

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

Returns the number of elements in row i.

Definition at line 91 of file sparsemat.cpp.

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

Definition at line 1819 of file sparsemat.cpp.

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

Definition at line 1760 of file sparsemat.cpp.

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

Definition at line 1788 of file sparsemat.cpp.

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

Definition at line 1506 of file sparsemat.cpp.

void mfem::SparseMatrix::SetDiagIdentity ( )

If a row contains only one diag entry of zero, set it to 1.

Definition at line 1138 of file sparsemat.cpp.

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

Definition at line 1700 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 1544 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 1577 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 160 of file sparsemat.cpp.

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

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

Definition at line 94 of file sparsemat.hpp.

void mfem::SparseMatrix::SortColumnIndices ( )

Sort the column indices corresponding to each row.

Definition at line 198 of file sparsemat.cpp.

void mfem::SparseMatrix::Symmetrize ( )

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

Definition at line 678 of file sparsemat.cpp.

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

Walks the sparse matrix.

Friends And Related Function Documentation

void Swap ( SparseMatrix A,
SparseMatrix B 
)
friend

Definition at line 2657 of file sparsemat.cpp.


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