12 #ifndef MFEM_SPARSEMAT_HPP
13 #define MFEM_SPARSEMAT_HPP
17 #include "../general/mem_alloc.hpp"
18 #include "../general/mem_manager.hpp"
19 #include "../general/device.hpp"
20 #include "../general/table.hpp"
21 #include "../general/globals.hpp"
28 #if defined(__alignas_is_defined)
72 #ifdef MFEM_USE_MEMALLOC
96 SparseMatrix(
int *i,
int *j,
double *data,
int m,
int n);
103 SparseMatrix(
int *i,
int *j,
double *data,
int m,
int n,
bool ownij,
104 bool owna,
bool issorted);
143 inline const int *
GetI()
const {
return I; }
148 inline const int *
GetJ()
const {
return J; }
156 int RowSize(
const int i)
const;
193 virtual double &
Elem(
int i,
int j);
196 virtual const double &
Elem(
int i,
int j)
const;
227 const double a = 1.0)
const;
256 const double a=1.0)
const;
356 double sc = 1.0)
const;
360 double sc = 1.0)
const;
370 void Finalize(
int skip_zeros,
bool fix_empty_rows);
378 void Threshold(
double tol,
bool fix_empty_rows =
false);
393 inline void SetColPtr(
const int row)
const;
404 inline void _Add_(
const int col,
const double a)
407 inline void _Set_(
const int col,
const double a)
410 inline double _Get_(
const int col)
const;
412 inline double &
SearchRow(
const int row,
const int col);
413 inline void _Add_(
const int row,
const int col,
const double a)
415 inline void _Set_(
const int row,
const int col,
const double a)
418 void Set(
const int i,
const int j,
const double a);
419 void Add(
const int i,
const int j,
const double a);
447 void ScaleRow(
const int row,
const double scale);
529 SparseMatrix *
Transpose(
const SparseMatrix &A);
540 SparseMatrix *
Mult(
const SparseMatrix &A,
const SparseMatrix &B,
541 SparseMatrix *OAB = NULL);
544 SparseMatrix *
TransposeMult(
const SparseMatrix &A,
const SparseMatrix &B);
548 const AbstractSparseMatrix &B);
551 DenseMatrix *
Mult(
const SparseMatrix &A, DenseMatrix &B);
554 DenseMatrix *
RAP(
const SparseMatrix &A, DenseMatrix &P);
557 DenseMatrix *
RAP(DenseMatrix &A,
const SparseMatrix &P);
561 SparseMatrix *
RAP(
const SparseMatrix &A,
const SparseMatrix &R,
562 SparseMatrix *ORAP = NULL);
565 SparseMatrix *
RAP(
const SparseMatrix &Rt,
const SparseMatrix &A,
566 const SparseMatrix &P);
569 SparseMatrix *
Mult_AtDA(
const SparseMatrix &A,
const Vector &D,
570 SparseMatrix *OAtDA = NULL);
574 SparseMatrix *
Add(
const SparseMatrix & A,
const SparseMatrix & B);
576 SparseMatrix *
Add(
double a,
const SparseMatrix & A,
double b,
577 const SparseMatrix & B);
579 SparseMatrix *
Add(Array<SparseMatrix *> & Ai);
582 void Add(
const SparseMatrix &A,
double alpha, DenseMatrix &B);
585 DenseMatrix *
OuterProduct(
const DenseMatrix &A,
const DenseMatrix &B);
588 SparseMatrix *
OuterProduct(
const DenseMatrix &A,
const SparseMatrix &B);
591 SparseMatrix *
OuterProduct(
const SparseMatrix &A,
const DenseMatrix &B);
594 SparseMatrix *
OuterProduct(
const SparseMatrix &A,
const SparseMatrix &B);
606 for (
int i = 0; i <
width; i++)
611 for (RowNode *node_p =
Rows[row]; node_p != NULL; node_p = node_p->Prev)
621 for (
int i = 0; i <
width; i++)
626 for (
int j =
I[row], end =
I[row+1]; j < end; j++)
639 node_p = node_p->Prev)
646 for (
int j =
I[
current_row], end =
I[current_row+1]; j < end; j++)
660 #ifdef MFEM_USE_MEMALLOC
663 node_p =
new RowNode;
666 node_p->Column = col;
670 return node_p->Value;
675 MFEM_VERIFY(j != -1,
"Entry for column " << col <<
" is not allocated.");
685 return (node_p == NULL) ? 0.0 : node_p->Value;
690 return (j == -1) ? 0.0 :
A[j];
700 for (node_p =
Rows[row]; 1; node_p = node_p->Prev)
704 #ifdef MFEM_USE_MEMALLOC
707 node_p =
new RowNode;
709 node_p->Prev =
Rows[row];
710 node_p->Column = col;
715 else if (node_p->Column == col)
720 return node_p->Value;
724 int *Ip =
I+row, *Jp =
J;
725 for (
int k = Ip[0], end = Ip[1]; k < end; k++)
732 MFEM_ABORT(
"Could not find entry for row = " << row <<
", col = " << col);
Memory< int > I
Array with size (height+1) containing the row offsets.
void _Add_(const int row, const int col, const double a)
double GetJacobiScaling() const
Determine appropriate scaling for Jacobi iteration.
int CheckFinite() const
Count the number of entries that are NOT finite, i.e. Inf or Nan.
static MemoryClass GetMemoryClass()
Get the current Device MemoryClass. This is the MemoryClass used by most MFEM device kernels to acces...
SparseMatrix * TransposeAbstractSparseMatrix(const AbstractSparseMatrix &A, int useActualWidth)
Transpose of a sparse matrix. A does not need to be a CSR matrix.
int RowSize(const int i) const
Returns the number of elements in row i.
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
void Add(const int i, const int j, const double a)
void _Add_(const int col, const double a)
Add a value to an entry in the "current row". See SetColPtr().
void BuildTranspose() const
Build and store internally the transpose of this matrix which will be used in the methods AddMultTran...
void Clear()
Clear the contents of the SparseMatrix.
void EliminateCol(int col, DiagonalPolicy dpolicy=DIAG_ZERO)
Eliminates the column col from the matrix.
SparseMatrix * At
Transpose of A. Owned. Used to perform MultTranspose() on devices.
void _Set_(const int row, const int col, const double a)
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
void DiagScale(const Vector &b, Vector &x, double sc=1.0) const
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
void MakeRef(const SparseMatrix &master)
Clear the contents of the SparseMatrix and make it a reference to master.
double & SearchRow(const int col)
Perform a fast search for an entry in the "current row". See SetColPtr().
DenseMatrix * ToDenseMatrix() const
Produces a DenseMatrix from a SparseMatrix.
bool Empty() const
Check if the SparseMatrix is empty.
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Set the diagonal value to zero.
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
int * GetJ()
Return the array J.
Abstract data type for sparse matrices.
SparseMatrix & operator+=(const SparseMatrix &B)
Add the sparse matrix 'B' to '*this'. This operation will cause an error if '*this' is finalized and ...
Data type dense matrix using column-major storage.
void AddMult(const Vector &x, Vector &y, const double a=1.0) const
y += A * x (default) or y += a * A * x
virtual ~SparseMatrix()
Destroys sparse matrix.
int * GetI()
Return the array I.
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
y += At * x (default) or y += a * At * x
Abstract data type for matrix inverse.
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 elim...
void PrintInfo(std::ostream &out) const
Print various sparse matrix staticstics.
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
void LoseData()
Lose the ownership of the graph (I, J) and data (A) arrays.
void AddRow(const int row, const Array< int > &cols, const Vector &srow)
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
void SortColumnIndices()
Sort the column indices corresponding to each row.
void PrintMM(std::ostream &out=mfem::out) const
Prints matrix in Matrix Market sparse format.
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const
void _Set_(const int col, const double a)
Set an entry in the "current row". See SetColPtr().
void BooleanMultTranspose(const Array< int > &x, Array< int > &y) const
y = At * x, treating all entries as booleans (zero=false, nonzero=true).
SparseMatrix * Mult_AtDA(const SparseMatrix &A, const Vector &D, SparseMatrix *OAtDA)
Matrix multiplication A^t D A. All matrices must be finalized.
void PartAddMult(const Array< int > &rows, const Vector &x, Vector &y, const double a=1.0) const
bool RowIsEmpty(const int row) const
void GetBlocks(Array2D< SparseMatrix * > &blocks) const
void ScaleRow(const int row, const double scale)
double * GetData()
Return the element data, i.e. the array A.
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
void MoveDiagonalFirst()
Move the diagonal entry to the first position in each row, preserving the order of the rest of the co...
void SetGraphOwner(bool ownij)
Set the graph ownership flag (I and J arrays).
ID for class SparseMatrix.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
bool OwnsData() const
Get the data ownership flag (A array).
void ClearColPtr() const
Reset the "current row" set by calling SetColPtr(). This method must be called between any two calls ...
void ScaleColumns(const Vector &sr)
this = this * diag(sr);
void PrintCSR2(std::ostream &out) const
Prints a sparse matrix to stream out in CSR format.
virtual int GetRow(const int row, Array< int > &cols, Vector &srow) const
Extract all column indices and values from a given row.
bool isSorted
Are the columns sorted already.
Memory< double > A
Array with size I[height], containing the actual entries of the sparse matrix, as indexed by the I ar...
SparseMatrix()
Create an empty SparseMatrix.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
bool OwnsGraph() const
Get the graph ownership flag (I and J arrays).
void SetWidth(int width_=-1)
Change the width of a SparseMatrix.
bool areColumnsSorted() const
double IsSymmetric() const
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix.
void SetDiagIdentity()
If a row contains only one diag entry of zero, set it to 1.
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
bool OwnsHostPtr() const
Return true if the host pointer is owned. Ownership indicates whether the pointer will be deleted by ...
void ScaleRows(const Vector &sl)
this = diag(sl) * this;
void Swap< SparseMatrix >(SparseMatrix &a, SparseMatrix &b)
Specialization of the template function Swap<> for class SparseMatrix.
DenseMatrix * OuterProduct(const DenseMatrix &A, const DenseMatrix &B)
Produces a block matrix with blocks A_{ij}*B.
void PrintMatlab(std::ostream &out=mfem::out) const
Prints matrix in matlab format.
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Set the diagonal value to one.
Dynamic 2D array using row-major layout.
void Jacobi(const Vector &b, const Vector &x0, Vector &x1, double sc) const
void SetSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
const double * GetData() const
Return the element data, i.e. the array A, const version.
Type
Enumeration defining IDs for some classes derived from Operator.
void SparseMatrixFunction(SparseMatrix &S, double(*f)(double))
Applies f() to each element of the matrix (after it is finalized).
virtual MatrixInverse * Inverse() const
This virtual method is not supported: it always returns NULL.
void EliminateRowColDiag(int rc, double value)
Perform elimination and set the diagonal entry to the given value.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void SetColPtr(const int row) const
Initialize the SparseMatrix for fast access to the entries of the given row which becomes the "curren...
const int * GetJ() const
Return the array J, const version.
void SetHostPtrOwner(bool own) const
Set/clear the ownership flag for the host pointer. Ownership indicates whether the pointer will be de...
void SetSubMatrixTranspose(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
void Set(const int i, const int j, const double a)
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.
void Threshold(double tol, bool fix_empty_rows=false)
Remove entries smaller in absolute value than a given tolerance tol. If fix_empty_rows is true...
int Size() const
For backward compatibility, define Size() to be synonym of Height().
void Gauss_Seidel_back(const Vector &x, Vector &y) const
void Gauss_Seidel_forw(const Vector &x, Vector &y) const
Gauss-Seidel forward and backward iterations over a vector x.
RowNode ** Rows
Array of linked lists, one for every row. This array represents the linked list (LIL) storage format...
void EliminateRow(int row, const double sol, Vector &rhs)
Eliminates a column from the transpose matrix.
int ActualWidth() const
Returns the actual Width of the matrix.
void GetDiag(Vector &d) const
Returns the Diagonal of A.
MemAlloc< RowNode, 1024 > RowNodeAlloc
Host memory; using new[] and delete[].
double GetRowNorml1(int irow) const
For i = irow compute .
void ResetTranspose() const
SparseMatrix * TransposeMult(const SparseMatrix &A, const SparseMatrix &B)
C = A^T B.
virtual void EliminateZeroRows(const double threshold=1e-12)
If a row contains only zeros, set its diagonal to 1.
bool Empty() const
Return true if the Memory object is empty, see Reset().
double & operator()(int i, int j)
Returns reference to A[i][j].
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
void EliminateCols(const Array< int > &cols, const Vector *x=NULL, Vector *b=NULL)
Eliminate all columns i for which cols[i] != 0.
SparseMatrix & operator=(const SparseMatrix &rhs)
Assignment operator: deep copy.
SparseMatrix & operator*=(double a)
void SetDataOwner(bool owna)
Set the data ownership flag (A array).
void SetRow(const int row, const Array< int > &cols, const Vector &srow)
void Print(std::ostream &out=mfem::out, int width_=4) const
Prints matrix to stream out.
const int * GetI() const
Return the array I, const version.
void Jacobi2(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
int MaxRowSize() const
Returns the maximum number of elements among all rows.
void Jacobi3(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
SparseMatrix * MultAbstractSparseMatrix(const AbstractSparseMatrix &A, const AbstractSparseMatrix &B)
Matrix product of sparse matrices. A and B do not need to be CSR matrices.
double _Get_(const int col) const
Read the value of an entry in the "current row". See SetColPtr().
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Memory< int > J
Array with size I[height], containing the column indices for all matrix entries, as indexed by the I ...
void Swap(SparseMatrix &other)
void PrintCSR(std::ostream &out) const
Prints matrix to stream out in hypre_CSRMatrix format.
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
MemoryClass
Memory classes identify subsets of memory types.
int width
Dimension of the input / number of columns in the matrix.
int CountSmallElems(double tol) const
Count the number of entries with |a_ij| <= tol.
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const
void GetRowSums(Vector &x) const
For all i compute .