MFEM  v3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
sparsemat.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_SPARSEMAT
13 #define MFEM_SPARSEMAT
14 
15 // Data types for sparse matrix
16 
17 #include "../general/mem_alloc.hpp"
18 #include "../general/table.hpp"
19 #include "densemat.hpp"
20 #include <iostream>
21 
22 namespace mfem
23 {
24 
25 class
26 #if defined(__alignas_is_defined)
27  alignas(double)
28 #endif
29  RowNode
30 {
31 public:
32  double Value;
33  RowNode *Prev;
34  int Column;
35 };
36 
39 {
40 private:
41 
45  int *I, *J;
46 
48  double *A;
49 
51  RowNode **Rows;
52 
53  mutable int current_row;
54  mutable int* ColPtrJ;
55  mutable RowNode ** ColPtrNode;
56 
57 #ifdef MFEM_USE_MEMALLOC
59  RowNodeAlloc * NodesMem;
60 #endif
61 
63  bool ownGraph;
65  bool ownData;
66 
68  bool isSorted;
69 
70  void Destroy(); // Delete all owned data
71  void SetEmpty(); // Init all entries with empty values
72 
73 public:
75  SparseMatrix() { SetEmpty(); }
76 
81  explicit SparseMatrix(int nrows, int ncols = -1);
82 
85  SparseMatrix(int *i, int *j, double *data, int m, int n);
86 
89  SparseMatrix(int *i, int *j, double *data, int m, int n, bool ownij,
90  bool owna, bool issorted);
91 
95  SparseMatrix(int nrows, int ncols, int rowsize);
96 
100  SparseMatrix(const SparseMatrix &mat, bool copy_graph = true);
101 
105  void MakeRef(const SparseMatrix &master);
106 
108  int Size() const { return Height(); }
109 
111  void Clear() { Destroy(); SetEmpty(); }
112 
114  bool Empty() { return (A == NULL) && (Rows == NULL); }
115 
117  inline int *GetI() const { return I; }
119  inline int *GetJ() const { return J; }
121  inline double *GetData() const { return A; }
123  int RowSize(const int i) const;
125  int MaxRowSize() const;
127  int *GetRowColumns(const int row);
128  const int *GetRowColumns(const int row) const;
130  double *GetRowEntries(const int row);
131  const double *GetRowEntries(const int row) const;
132 
134 
141  void SetWidth(int width_ = -1);
142 
144 
145  int ActualWidth();
146 
148  void SortColumnIndices();
149 
151  virtual double &Elem(int i, int j);
152 
154  virtual const double &Elem(int i, int j) const;
155 
157  double &operator()(int i, int j);
158 
160  const double &operator()(int i, int j) const;
161 
163  void GetDiag(Vector & d) const;
164 
166  virtual void Mult(const Vector &x, Vector &y) const;
167 
169  void AddMult(const Vector &x, Vector &y, const double a = 1.0) const;
170 
172  void MultTranspose(const Vector &x, Vector &y) const;
173 
175  void AddMultTranspose(const Vector &x, Vector &y,
176  const double a = 1.0) const;
177 
178  void PartMult(const Array<int> &rows, const Vector &x, Vector &y) const;
179  void PartAddMult(const Array<int> &rows, const Vector &x, Vector &y,
180  const double a=1.0) const;
181 
183  void BooleanMult(const Array<int> &x, Array<int> &y) const;
185  void BooleanMultTranspose(const Array<int> &x, Array<int> &y) const;
186 
188  double InnerProduct(const Vector &x, const Vector &y) const;
189 
191  void GetRowSums(Vector &x) const;
193  double GetRowNorml1(int irow) const;
194 
196  virtual MatrixInverse *Inverse() const;
197 
199  void EliminateRow(int row, const double sol, Vector &rhs);
201 
207  void EliminateRow(int row, int setOneDiagonal = 0);
208  void EliminateCol(int col);
210  void EliminateCols(Array<int> &cols, Vector *x = NULL, Vector *b = NULL);
211 
216  void EliminateRowCol(int rc, const double sol, Vector &rhs, int d = 0);
217 
220  void EliminateRowColMultipleRHS(int rc, const Vector &sol,
221  DenseMatrix &rhs, int d = 0);
222 
223  void EliminateRowCol(int rc, int d = 0);
225  void EliminateRowColDiag(int rc, double value);
226  // Same as above + save the eliminated entries in Ae so that
227  // (*this) + Ae is the original matrix
228  void EliminateRowCol(int rc, SparseMatrix &Ae, int d = 0);
229 
231  void SetDiagIdentity();
233  void EliminateZeroRows();
234 
236  void Gauss_Seidel_forw(const Vector &x, Vector &y) const;
237  void Gauss_Seidel_back(const Vector &x, Vector &y) const;
238 
240  double GetJacobiScaling() const;
243  void Jacobi(const Vector &b, const Vector &x0, Vector &x1, double sc) const;
244 
245  void DiagScale(const Vector &b, Vector &x, double sc = 1.0) const;
246 
248  void Jacobi2(const Vector &b, const Vector &x0, Vector &x1,
249  double sc = 1.0) const;
250 
252  void Jacobi3(const Vector &b, const Vector &x0, Vector &x1,
253  double sc = 1.0) const;
254 
259  virtual void Finalize(int skip_zeros = 1);
260 
261  bool Finalized() const { return (A != NULL); }
262  bool areColumnsSorted() const { return isSorted; }
263 
267  void GetBlocks(Array2D<SparseMatrix *> &blocks) const;
268 
269  void GetSubMatrix(const Array<int> &rows, const Array<int> &cols,
270  DenseMatrix &subm);
271 
272  inline void SetColPtr(const int row) const;
273  inline void ClearColPtr() const;
274  inline double &SearchRow(const int col);
275  inline void _Add_(const int col, const double a)
276  { SearchRow(col) += a; }
277  inline void _Set_(const int col, const double a)
278  { SearchRow(col) = a; }
279  inline double _Get_(const int col) const;
280 
281  inline double &SearchRow(const int row, const int col);
282  inline void _Add_(const int row, const int col, const double a)
283  { SearchRow(row, col) += a; }
284  inline void _Set_(const int row, const int col, const double a)
285  { SearchRow(row, col) = a; }
286 
287  void Set(const int i, const int j, const double a);
288  void Add(const int i, const int j, const double a);
289 
290  void SetSubMatrix(const Array<int> &rows, const Array<int> &cols,
291  const DenseMatrix &subm, int skip_zeros = 1);
292 
293  void SetSubMatrixTranspose(const Array<int> &rows, const Array<int> &cols,
294  const DenseMatrix &subm, int skip_zeros = 1);
295 
296  void AddSubMatrix(const Array<int> &rows, const Array<int> &cols,
297  const DenseMatrix &subm, int skip_zeros = 1);
298 
299  bool RowIsEmpty(const int row) const;
300 
309  virtual int GetRow(const int row, Array<int> &cols, Vector &srow) const;
310 
311  void SetRow(const int row, const Array<int> &cols, const Vector &srow);
312  void AddRow(const int row, const Array<int> &cols, const Vector &srow);
313 
314  void ScaleRow(const int row, const double scale);
315  // this = diag(sl) * this;
316  void ScaleRows(const Vector & sl);
317  // this = this * diag(sr);
318  void ScaleColumns(const Vector & sr);
319 
323 
326  void Add(const double a, const SparseMatrix &B);
327 
328  SparseMatrix &operator=(double a);
329 
330  SparseMatrix &operator*=(double a);
331 
333  void Print(std::ostream &out = std::cout, int width_ = 4) const;
334 
336  void PrintMatlab(std::ostream &out = std::cout) const;
337 
339  void PrintMM(std::ostream &out = std::cout) const;
340 
342  void PrintCSR(std::ostream &out) const;
343 
345  void PrintCSR2(std::ostream &out) const;
346 
348  int Walk(int &i, int &j, double &a);
349 
351  double IsSymmetric() const;
352 
354  void Symmetrize();
355 
357  virtual int NumNonZeroElems() const;
358 
359  double MaxNorm() const;
360 
362  int CountSmallElems(double tol) const;
363 
365  void SetGraphOwner(bool ownij) { ownGraph = ownij; }
367  void SetDataOwner(bool owna) { ownData = owna; }
369  bool OwnsGraph() const { return ownGraph; }
371  bool OwnsData() const { return ownData; }
373  void LoseData() { ownGraph = ownData = false; }
374 
375  void Swap(SparseMatrix &other);
376 
378  virtual ~SparseMatrix() { Destroy(); }
379 };
380 
382 void SparseMatrixFunction(SparseMatrix &S, double (*f)(double));
383 
384 
386 SparseMatrix *Transpose(const SparseMatrix &A);
388 SparseMatrix *TransposeAbstractSparseMatrix (const AbstractSparseMatrix &A,
389  int useActualWidth);
390 
397 SparseMatrix *Mult(const SparseMatrix &A, const SparseMatrix &B,
398  SparseMatrix *OAB = NULL);
399 
401 SparseMatrix *MultAbstractSparseMatrix (const AbstractSparseMatrix &A,
402  const AbstractSparseMatrix &B);
403 
405 DenseMatrix *Mult(const SparseMatrix &A, DenseMatrix &B);
406 
408 DenseMatrix *RAP(const SparseMatrix &A, DenseMatrix &P);
409 
412 SparseMatrix *RAP(const SparseMatrix &A, const SparseMatrix &R,
413  SparseMatrix *ORAP = NULL);
414 
416 SparseMatrix *RAP(const SparseMatrix &Rt, const SparseMatrix &A,
417  const SparseMatrix &P);
418 
420 SparseMatrix *Mult_AtDA(const SparseMatrix &A, const Vector &D,
421  SparseMatrix *OAtDA = NULL);
422 
423 
425 SparseMatrix * Add(const SparseMatrix & A, const SparseMatrix & B);
427 SparseMatrix * Add(double a, const SparseMatrix & A, double b,
428  const SparseMatrix & B);
430 SparseMatrix * Add(Array<SparseMatrix *> & Ai);
431 
432 
433 // Inline methods
434 
435 inline void SparseMatrix::SetColPtr(const int row) const
436 {
437  if (Rows)
438  {
439  if (ColPtrNode == NULL)
440  {
441  ColPtrNode = new RowNode *[width];
442  for (int i = 0; i < width; i++)
443  {
444  ColPtrNode[i] = NULL;
445  }
446  }
447  for (RowNode *node_p = Rows[row]; node_p != NULL; node_p = node_p->Prev)
448  {
449  ColPtrNode[node_p->Column] = node_p;
450  }
451  }
452  else
453  {
454  if (ColPtrJ == NULL)
455  {
456  ColPtrJ = new int[width];
457  for (int i = 0; i < width; i++)
458  {
459  ColPtrJ[i] = -1;
460  }
461  }
462  for (int j = I[row], end = I[row+1]; j < end; j++)
463  {
464  ColPtrJ[J[j]] = j;
465  }
466  }
467  current_row = row;
468 }
469 
470 inline void SparseMatrix::ClearColPtr() const
471 {
472  if (Rows)
473  for (RowNode *node_p = Rows[current_row]; node_p != NULL;
474  node_p = node_p->Prev)
475  {
476  ColPtrNode[node_p->Column] = NULL;
477  }
478  else
479  for (int j = I[current_row], end = I[current_row+1]; j < end; j++)
480  {
481  ColPtrJ[J[j]] = -1;
482  }
483 }
484 
485 inline double &SparseMatrix::SearchRow(const int col)
486 {
487  if (Rows)
488  {
489  RowNode *node_p = ColPtrNode[col];
490  if (node_p == NULL)
491  {
492 #ifdef MFEM_USE_MEMALLOC
493  node_p = NodesMem->Alloc();
494 #else
495  node_p = new RowNode;
496 #endif
497  node_p->Prev = Rows[current_row];
498  node_p->Column = col;
499  node_p->Value = 0.0;
500  Rows[current_row] = ColPtrNode[col] = node_p;
501  }
502  return node_p->Value;
503  }
504  else
505  {
506  const int j = ColPtrJ[col];
507  MFEM_VERIFY(j != -1, "Entry for column " << col << " is not allocated.");
508  return A[j];
509  }
510 }
511 
512 inline double SparseMatrix::_Get_(const int col) const
513 {
514  if (Rows)
515  {
516  RowNode *node_p = ColPtrNode[col];
517  return (node_p == NULL) ? 0.0 : node_p->Value;
518  }
519  else
520  {
521  const int j = ColPtrJ[col];
522  return (j == -1) ? 0.0 : A[j];
523  }
524 }
525 
526 inline double &SparseMatrix::SearchRow(const int row, const int col)
527 {
528  if (Rows)
529  {
530  RowNode *node_p;
531 
532  for (node_p = Rows[row]; 1; node_p = node_p->Prev)
533  {
534  if (node_p == NULL)
535  {
536 #ifdef MFEM_USE_MEMALLOC
537  node_p = NodesMem->Alloc();
538 #else
539  node_p = new RowNode;
540 #endif
541  node_p->Prev = Rows[row];
542  node_p->Column = col;
543  node_p->Value = 0.0;
544  Rows[row] = node_p;
545  break;
546  }
547  else if (node_p->Column == col)
548  {
549  break;
550  }
551  }
552  return node_p->Value;
553  }
554  else
555  {
556  int *Ip = I+row, *Jp = J;
557  for (int k = Ip[0], end = Ip[1]; k < end; k++)
558  {
559  if (Jp[k] == col)
560  {
561  return A[k];
562  }
563  }
564  MFEM_ABORT("Could not find entry for row = " << row << ", col = " << col);
565  }
566  return A[0];
567 }
568 
570 template<> inline void Swap<SparseMatrix>(SparseMatrix &a, SparseMatrix &b)
571 {
572  a.Swap(b);
573 }
574 
575 }
576 
577 #endif
Elem * Alloc()
Definition: mem_alloc.hpp:124
void _Add_(const int row, const int col, const double a)
Definition: sparsemat.hpp:282
double GetJacobiScaling() const
Determine appropriate scaling for Jacobi iteration.
Definition: sparsemat.cpp:1615
SparseMatrix * TransposeAbstractSparseMatrix(const AbstractSparseMatrix &A, int useActualWidth)
Transpose of a sparse matrix. A does not need to be a CSR matrix.
Definition: sparsemat.cpp:2517
int RowSize(const int i) const
Returns the number of elements in row i.
Definition: sparsemat.cpp:212
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:905
void Add(const int i, const int j, const double a)
Definition: sparsemat.cpp:1816
void _Add_(const int col, const double a)
Definition: sparsemat.hpp:275
void Clear()
Clear the contents of the SparseMatrix.
Definition: sparsemat.hpp:111
void _Set_(const int row, const int col, const double a)
Definition: sparsemat.hpp:284
void Print(std::ostream &out=std::cout, int width_=4) const
Prints matrix to stream out.
Definition: sparsemat.cpp:2247
void DiagScale(const Vector &b, Vector &x, double sc=1.0) const
Definition: sparsemat.cpp:1679
void MakeRef(const SparseMatrix &master)
Definition: sparsemat.cpp:186
double & SearchRow(const int col)
Definition: sparsemat.hpp:485
void PrintMM(std::ostream &out=std::cout) const
Prints matrix in Matrix Market sparse format.
Definition: sparsemat.cpp:2309
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Definition: sparsemat.cpp:259
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, but treat all elements as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:595
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:451
Abstract data type for sparse matrices.
Definition: matrix.hpp:69
SparseMatrix & operator+=(const SparseMatrix &B)
Definition: sparsemat.cpp:2158
void EliminateRowCol(int rc, const double sol, Vector &rhs, int d=0)
Definition: sparsemat.cpp:1071
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
void AddMult(const Vector &x, Vector &y, const double a=1.0) const
y += A * x (default) or y += a * A * x
Definition: sparsemat.cpp:445
double MaxNorm() const
Definition: sparsemat.cpp:925
virtual ~SparseMatrix()
Destroys sparse matrix.
Definition: sparsemat.hpp:378
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const
y += At * x (default) or y += a * At * x
Definition: sparsemat.cpp:520
Abstract data type for matrix inverse.
Definition: matrix.hpp:58
void EliminateRowColMultipleRHS(int rc, const Vector &sol, DenseMatrix &rhs, int d=0)
Definition: sparsemat.cpp:1151
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
Definition: sparsemat.cpp:354
void LoseData()
Lose the ownership of the graph (I, J) and data (A) arrays.
Definition: sparsemat.hpp:373
void AddRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:2039
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
Definition: sparsemat.cpp:273
void SortColumnIndices()
Sort the column indices corresponding to each row.
Definition: sparsemat.cpp:325
void _Set_(const int col, const double a)
Definition: sparsemat.hpp:277
void BooleanMultTranspose(const Array< int > &x, Array< int > &y) const
y = At * x, but treat all elements as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:618
SparseMatrix * Mult_AtDA(const SparseMatrix &A, const Vector &D, SparseMatrix *OAtDA)
Matrix multiplication A^t D A. All matrices must be finalized.
Definition: sparsemat.cpp:2862
void PartAddMult(const Array< int > &rows, const Vector &x, Vector &y, const double a=1.0) const
Definition: sparsemat.cpp:577
bool RowIsEmpty(const int row) const
Definition: sparsemat.cpp:1931
void GetBlocks(Array2D< SparseMatrix * > &blocks) const
Definition: sparsemat.cpp:778
void ScaleRow(const int row, const double scale)
Definition: sparsemat.cpp:2071
int ActualWidth()
Returns the actual Width of the matrix.
Definition: sparsemat.cpp:2426
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
Definition: sparsemat.cpp:890
void SetGraphOwner(bool ownij)
Set the graph ownership flag (I and J arrays).
Definition: sparsemat.hpp:365
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2809
bool OwnsData() const
Get the data ownership flag (A array).
Definition: sparsemat.hpp:371
void ClearColPtr() const
Definition: sparsemat.hpp:470
void EliminateCol(int col)
Definition: sparsemat.cpp:1026
void ScaleColumns(const Vector &sr)
Definition: sparsemat.cpp:2130
void PrintCSR2(std::ostream &out) const
Prints a sparse matrix to stream out in CSR format.
Definition: sparsemat.cpp:2353
virtual int GetRow(const int row, Array< int > &cols, Vector &srow) const
Definition: sparsemat.cpp:1952
SparseMatrix()
Create an empty SparseMatrix.
Definition: sparsemat.hpp:75
Data type sparse matrix.
Definition: sparsemat.hpp:38
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
Definition: operator.hpp:35
void EliminateCols(Array< int > &cols, Vector *x=NULL, Vector *b=NULL)
Eliminate all columns &#39;i&#39; for which cols[i] != 0.
Definition: sparsemat.cpp:1040
bool OwnsGraph() const
Get the graph ownership flag (I and J arrays).
Definition: sparsemat.hpp:369
void SetWidth(int width_=-1)
Change the width of a SparseMatrix.
Definition: sparsemat.cpp:287
bool areColumnsSorted() const
Definition: sparsemat.hpp:262
HypreParMatrix * RAP(HypreParMatrix *A, HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1365
double IsSymmetric() const
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix.
Definition: sparsemat.cpp:868
void SetDiagIdentity()
If a row contains only one diag entry of zero, set it to 1.
Definition: sparsemat.cpp:1429
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
Definition: table.cpp:391
void ScaleRows(const Vector &sl)
Definition: sparsemat.cpp:2099
void Swap< SparseMatrix >(SparseMatrix &a, SparseMatrix &b)
Specialization of the template function Swap&lt;&gt; for class SparseMatrix.
Definition: sparsemat.hpp:570
virtual void Finalize(int skip_zeros=1)
Definition: sparsemat.cpp:706
void EliminateZeroRows()
If a row contains only zeros, set its diagonal to 1.
Definition: sparsemat.cpp:1438
double * GetData() const
Return element data.
Definition: sparsemat.hpp:121
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:1759
int * GetI() const
Return the array I.
Definition: sparsemat.hpp:117
void Jacobi(const Vector &b, const Vector &x0, Vector &x1, double sc) const
Definition: sparsemat.cpp:1648
void SetSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:1835
bool Finalized() const
Definition: sparsemat.hpp:261
void SparseMatrixFunction(SparseMatrix &S, double(*f)(double))
Applies f() to each element of the matrix (after it is finalized).
Definition: sparsemat.cpp:2452
virtual MatrixInverse * Inverse() const
Returns a pointer to approximation of the matrix inverse.
Definition: sparsemat.cpp:978
void EliminateRowColDiag(int rc, double value)
Perform elimination and set the diagonal entry to the given value.
Definition: sparsemat.cpp:1301
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:439
void SetColPtr(const int row) const
Definition: sparsemat.hpp:435
SparseMatrix & operator=(double a)
Definition: sparsemat.cpp:2211
void SetSubMatrixTranspose(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:1868
void Set(const int i, const int j, const double a)
Definition: sparsemat.cpp:1797
void PrintMatlab(std::ostream &out=std::cout) const
Prints matrix in matlab format.
Definition: sparsemat.cpp:2291
int Size() const
For backward compatibility define Size to be synonym of Height()
Definition: sparsemat.hpp:108
void Gauss_Seidel_back(const Vector &x, Vector &y) const
Definition: sparsemat.cpp:1540
void Gauss_Seidel_forw(const Vector &x, Vector &y) const
Gauss-Seidel forward and backward iterations over a vector x.
Definition: sparsemat.cpp:1465
void EliminateRow(int row, const double sol, Vector &rhs)
Eliminates a column from the transpose matrix.
Definition: sparsemat.cpp:983
void GetDiag(Vector &d) const
Returns the Diagonal of A.
Definition: sparsemat.cpp:411
double GetRowNorml1(int irow) const
For i = irow compute .
Definition: sparsemat.cpp:686
double & operator()(int i, int j)
Returns reference to A[i][j].
Definition: sparsemat.cpp:364
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:514
SparseMatrix & operator*=(double a)
Definition: sparsemat.cpp:2229
void SetDataOwner(bool owna)
Set the data ownership flag (A array).
Definition: sparsemat.hpp:367
bool Empty()
Check if the SparseMatrix is empty.
Definition: sparsemat.hpp:114
void SetRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:1991
Vector data type.
Definition: vector.hpp:33
int Walk(int &i, int &j, double &a)
Walks the sparse matrix.
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm)
Definition: sparsemat.cpp:1903
void Jacobi2(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
Definition: sparsemat.cpp:1711
int MaxRowSize() const
Returns the maximum number of elements among all rows.
Definition: sparsemat.cpp:235
void Jacobi3(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
Definition: sparsemat.cpp:1735
SparseMatrix * MultAbstractSparseMatrix(const AbstractSparseMatrix &A, const AbstractSparseMatrix &B)
Matrix product of sparse matrices. A and B do not need to be CSR matrices.
Definition: sparsemat.cpp:2720
double _Get_(const int col) const
Definition: sparsemat.hpp:512
void Swap(SparseMatrix &other)
Definition: sparsemat.cpp:2990
int * GetJ() const
Return the array J.
Definition: sparsemat.hpp:119
void PrintCSR(std::ostream &out) const
Prints matrix to stream out in hypre_CSRMatrix format.
Definition: sparsemat.cpp:2329
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:641
int CountSmallElems(double tol) const
Count the number of entries with |a_ij| &lt; tol.
Definition: sparsemat.cpp:948
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const
Definition: sparsemat.cpp:559
void GetRowSums(Vector &x) const
For all i compute .
Definition: sparsemat.cpp:667