12 #ifndef MFEM_SPARSEMAT
13 #define MFEM_SPARSEMAT
17 #include "../general/mem_alloc.hpp"
18 #include "../general/table.hpp"
26 #if defined(__alignas_is_defined)
53 mutable int current_row;
55 mutable RowNode ** ColPtrNode;
57 #ifdef MFEM_USE_MEMALLOC
85 SparseMatrix(
int *i,
int *j,
double *data,
int m,
int n);
89 SparseMatrix(
int *i,
int *j,
double *data,
int m,
int n,
bool ownij,
90 bool owna,
bool issorted);
106 void Clear() { Destroy(); SetEmpty(); }
109 bool Empty() {
return (A == NULL) && (Rows == NULL); }
112 inline int *
GetI()
const {
return I; }
114 inline int *
GetJ()
const {
return J; }
118 int RowSize(
const int i)
const;
146 virtual double &
Elem(
int i,
int j);
149 virtual const double &
Elem(
int i,
int j)
const;
171 const double a = 1.0)
const;
175 const double a=1.0)
const;
244 double sc = 1.0)
const;
248 double sc = 1.0)
const;
254 virtual void Finalize(
int skip_zeros = 1);
267 inline void SetColPtr(
const int row)
const;
270 inline void _Add_(
const int col,
const double a)
272 inline void _Set_(
const int col,
const double a)
274 inline double _Get_(
const int col)
const;
276 inline double &
SearchRow(
const int row,
const int col);
277 inline void _Add_(
const int row,
const int col,
const double a)
279 inline void _Set_(
const int row,
const int col,
const double a)
282 void Set(
const int i,
const int j,
const double a);
283 void Add(
const int i,
const int j,
const double a);
309 void ScaleRow(
const int row,
const double scale);
328 void Print(std::ostream &out = std::cout,
int width_ = 4)
const;
331 void PrintMatlab(std::ostream &out = std::cout)
const;
334 void PrintMM(std::ostream &out = std::cout)
const;
337 void PrintCSR(std::ostream &out)
const;
343 int Walk(
int &i,
int &j,
double &a);
381 SparseMatrix *
Transpose(
const SparseMatrix &A);
392 SparseMatrix *
Mult(
const SparseMatrix &A,
const SparseMatrix &B,
393 SparseMatrix *OAB = NULL);
397 const AbstractSparseMatrix &B);
402 SparseMatrix *
RAP(
const SparseMatrix &A,
const SparseMatrix &R,
403 SparseMatrix *ORAP = NULL);
406 SparseMatrix *
RAP(
const SparseMatrix &Rt,
const SparseMatrix &A,
407 const SparseMatrix &P);
410 SparseMatrix *
Mult_AtDA(
const SparseMatrix &A,
const Vector &D,
411 SparseMatrix *OAtDA = NULL);
415 SparseMatrix *
Add(
const SparseMatrix & A,
const SparseMatrix & B);
417 SparseMatrix *
Add(
double a,
const SparseMatrix & A,
double b,
418 const SparseMatrix & B);
420 SparseMatrix *
Add(Array<SparseMatrix *> & Ai);
429 if (ColPtrNode == NULL)
431 ColPtrNode =
new RowNode *[
width];
432 for (
int i = 0; i <
width; i++)
434 ColPtrNode[i] = NULL;
437 for (RowNode *node_p = Rows[row]; node_p != NULL; node_p = node_p->Prev)
439 ColPtrNode[node_p->Column] = node_p;
446 ColPtrJ =
new int[
width];
447 for (
int i = 0; i <
width; i++)
452 for (
int j = I[row], end = I[row+1]; j < end; j++)
463 for (RowNode *node_p = Rows[current_row]; node_p != NULL;
464 node_p = node_p->Prev)
466 ColPtrNode[node_p->Column] = NULL;
469 for (
int j = I[current_row], end = I[current_row+1]; j < end; j++)
479 RowNode *node_p = ColPtrNode[col];
482 #ifdef MFEM_USE_MEMALLOC
483 node_p = NodesMem->
Alloc();
485 node_p =
new RowNode;
487 node_p->Prev = Rows[current_row];
488 node_p->Column = col;
490 Rows[current_row] = ColPtrNode[col] = node_p;
492 return node_p->Value;
496 const int j = ColPtrJ[col];
497 MFEM_VERIFY(j != -1,
"Entry for column " << col <<
" is not allocated.");
506 RowNode *node_p = ColPtrNode[col];
507 return (node_p == NULL) ? 0.0 : node_p->Value;
511 const int j = ColPtrJ[col];
512 return (j == -1) ? 0.0 : A[j];
522 for (node_p = Rows[row]; 1; node_p = node_p->Prev)
526 #ifdef MFEM_USE_MEMALLOC
527 node_p = NodesMem->
Alloc();
529 node_p =
new RowNode;
531 node_p->Prev = Rows[row];
532 node_p->Column = col;
537 else if (node_p->Column == col)
542 return node_p->Value;
546 int *Ip = I+row, *Jp = J;
547 for (
int k = Ip[0], end = Ip[1]; k < end; k++)
554 MFEM_ABORT(
"Could not find entry for row = " << row <<
", col = " << col);
void _Add_(const int row, const int col, const double a)
double GetJacobiScaling() const
Determine appropriate scaling for Jacobi iteration.
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)
void Clear()
Clear the contents of the SparseMatrix.
void _Set_(const int row, const int col, const double a)
void Print(std::ostream &out=std::cout, int width_=4) const
Prints matrix to stream out.
void DiagScale(const Vector &b, Vector &x, double sc=1.0) const
void MakeRef(const SparseMatrix &master)
double & SearchRow(const int col)
void PrintMM(std::ostream &out=std::cout) const
Prints matrix in Matrix Market sparse format.
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, but treat all elements as booleans (zero=false, nonzero=true).
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Abstract data type for sparse matrices.
void EliminateRowCol(int rc, const double sol, Vector &rhs, int d=0)
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.
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, int d=0)
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 _Set_(const int col, const double a)
void BooleanMultTranspose(const Array< int > &x, Array< int > &y) const
y = At * x, but treat all elements 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)
int ActualWidth()
Returns the actual Width of the matrix.
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
void SetGraphOwner(bool ownij)
Set the graph ownership flag (I and J arrays).
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 EliminateCol(int col)
void ScaleColumns(const Vector &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
SparseMatrix()
Create an empty SparseMatrix.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
void EliminateCols(Array< int > &cols, Vector *x=NULL, Vector *b=NULL)
Eliminate all columns 'i' for which cols[i] != 0.
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
HypreParMatrix * RAP(HypreParMatrix *A, HypreParMatrix *P)
Returns the matrix P^t * A * P.
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.
void ScaleRows(const Vector &sl)
void Swap< SparseMatrix >(SparseMatrix &a, SparseMatrix &b)
Specialization of the template function Swap<> for class SparseMatrix.
virtual void Finalize(int skip_zeros=1)
void EliminateZeroRows()
If a row contains only zeros, set its diagonal to 1.
double * GetData() const
Return element data.
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
int * GetI() const
Return the array I.
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)
void SparseMatrixFunction(SparseMatrix &S, double(*f)(double))
Applies f() to each element of the matrix (after it is finalized).
virtual MatrixInverse * Inverse() const
Returns a pointer to approximation of the matrix inverse.
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
SparseMatrix & operator=(double a)
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 PrintMatlab(std::ostream &out=std::cout) const
Prints matrix in matlab format.
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.
void EliminateRow(int row, const double sol, Vector &rhs)
Eliminates a column from the transpose matrix.
void GetDiag(Vector &d) const
Returns the Diagonal of A.
double GetRowNorml1(int irow) const
For i = irow compute .
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.
SparseMatrix & operator*=(double a)
void SetDataOwner(bool owna)
Set the data ownership flag (A array).
bool Empty()
Check if the SparseMatrix is empty.
void SetRow(const int row, const Array< int > &cols, const Vector &srow)
int Walk(int &i, int &j, double &a)
Walks the sparse matrix.
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm)
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
SparseMatrix & operator+=(SparseMatrix &B)
void Swap(SparseMatrix &other)
int * GetJ() const
Return the array J.
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.
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 .