MFEM  v4.2.0
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-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_SPARSEMAT_HPP
13 #define MFEM_SPARSEMAT_HPP
14 
15 // Data types for sparse matrix
16 
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"
22 #include "densemat.hpp"
23 
24 #ifdef MFEM_USE_CUDA
25 #include <cusparse.h>
26 #include <library_types.h>
27 #include "../general/cuda.hpp"
28 #endif
29 
30 namespace mfem
31 {
32 
33 class
34 #if defined(__alignas_is_defined)
35  alignas(double)
36 #endif
37  RowNode
38 {
39 public:
40  double Value;
41  RowNode *Prev;
42  int Column;
43 };
44 
45 /// Data type sparse matrix
47 {
48 protected:
49  /// @name Arrays used by the CSR storage format.
50  /** */
51  ///@{
52  /// @brief %Array with size (#height+1) containing the row offsets.
53  /** The data for row r, 0 <= r < height, is at offsets j, I[r] <= j < I[r+1].
54  The offsets, j, are indices in the #J and #A arrays. The first entry in
55  this array is always zero, I[0] = 0, and the last entry, I[height], gives
56  the total number of entries stored (at a minimum, all nonzeros must be
57  represented) in the sparse matrix. */
59  /** @brief %Array with size #I[#height], containing the column indices for
60  all matrix entries, as indexed by the #I array. */
62  /** @brief %Array with size #I[#height], containing the actual entries of the
63  sparse matrix, as indexed by the #I array. */
65  ///@}
66 
67  /** @brief %Array of linked lists, one for every row. This array represents
68  the linked list (LIL) storage format. */
69  RowNode **Rows;
70 
71  mutable int current_row;
72  mutable int* ColPtrJ;
73  mutable RowNode ** ColPtrNode;
74 
75  /// Transpose of A. Owned. Used to perform MultTranspose() on devices.
76  mutable SparseMatrix *At;
77 
78 #ifdef MFEM_USE_MEMALLOC
81 #endif
82 
83  /// Are the columns sorted already.
84  bool isSorted;
85 
86  void Destroy(); // Delete all owned data
87  void SetEmpty(); // Init all entries with empty values
88 
89  bool useCuSparse{true}; // Use cuSPARSE if available
90 
91  // Initialize cuSPARSE
92  void InitCuSparse();
93 
94 #ifdef MFEM_USE_CUDA
95  cusparseStatus_t status;
96  static cusparseHandle_t handle;
97  cusparseMatDescr_t descr=0;
98  static size_t bufferSize;
99  static void *dBuffer;
100  mutable bool initBuffers{false};
101 
102  static int SparseMatrixCount;
103  mutable cusparseSpMatDescr_t matA_descr;
104  mutable cusparseDnVecDescr_t vecX_descr;
105  mutable cusparseDnVecDescr_t vecY_descr;
106 #endif
107 
108 public:
109  /// Create an empty SparseMatrix.
111  {
112  SetEmpty();
113 
114  InitCuSparse();
115  }
116 
117  /** @brief Create a sparse matrix with flexible sparsity structure using a
118  row-wise linked list (LIL) format. */
119  /** New entries are added as needed by methods like AddSubMatrix(),
120  SetSubMatrix(), etc. Calling Finalize() will convert the SparseMatrix to
121  the more compact compressed sparse row (CSR) format. */
122  explicit SparseMatrix(int nrows, int ncols = -1);
123 
124  /** @brief Create a sparse matrix in CSR format. Ownership of @a i, @a j, and
125  @a data is transferred to the SparseMatrix. */
126  SparseMatrix(int *i, int *j, double *data, int m, int n);
127 
128  /** @brief Create a sparse matrix in CSR format. Ownership of @a i, @a j, and
129  @a data is optionally transferred to the SparseMatrix. */
130  /** If the parameter @a data is NULL, then the internal #A array is allocated
131  by this constructor (initializing it with zeros and taking ownership,
132  regardless of the parameter @a owna). */
133  SparseMatrix(int *i, int *j, double *data, int m, int n, bool ownij,
134  bool owna, bool issorted);
135 
136  /** @brief Create a sparse matrix in CSR format where each row has space
137  allocated for exactly @a rowsize entries. */
138  /** SetRow() can then be called or the #I, #J, #A arrays can be used
139  directly. */
140  SparseMatrix(int nrows, int ncols, int rowsize);
141 
142  /// Copy constructor (deep copy).
143  /** If @a mat is finalized and @a copy_graph is false, the #I and #J arrays
144  will use a shallow copy (copy the pointers only) without transferring
145  ownership. */
146  SparseMatrix(const SparseMatrix &mat, bool copy_graph = true);
147 
148  /// Create a SparseMatrix with diagonal @a v, i.e. A = Diag(v)
149  SparseMatrix(const Vector & v);
150 
151  // Runtime option to use cuSPARSE. Only valid when using a CUDA backend.
152  void UseCuSparse(bool _useCuSparse = true) { useCuSparse = _useCuSparse;}
153 
154  /// Assignment operator: deep copy
155  SparseMatrix& operator=(const SparseMatrix &rhs);
156 
157  /** @brief Clear the contents of the SparseMatrix and make it a reference to
158  @a master */
159  /** After this call, the matrix will point to the same data as @a master but
160  it will not own its data. The @a master must be finalized. */
161  void MakeRef(const SparseMatrix &master);
162 
163  /// For backward compatibility, define Size() to be synonym of Height().
164  int Size() const { return Height(); }
165 
166  /// Clear the contents of the SparseMatrix.
167  void Clear() { Destroy(); SetEmpty(); }
168 
169  /// Check if the SparseMatrix is empty.
170  bool Empty() const { return (A == NULL) && (Rows == NULL); }
171 
172  /// Return the array #I.
173  inline int *GetI() { return I; }
174  /// Return the array #I, const version.
175  inline const int *GetI() const { return I; }
176 
177  /// Return the array #J.
178  inline int *GetJ() { return J; }
179  /// Return the array #J, const version.
180  inline const int *GetJ() const { return J; }
181 
182  /// Return the element data, i.e. the array #A.
183  inline double *GetData() { return A; }
184  /// Return the element data, i.e. the array #A, const version.
185  inline const double *GetData() const { return A; }
186 
187  // Memory access methods for the #I array.
188  Memory<int> &GetMemoryI() { return I; }
189  const Memory<int> &GetMemoryI() const { return I; }
190  const int *ReadI(bool on_dev = true) const
191  { return mfem::Read(I, Height()+1, on_dev); }
192  int *WriteI(bool on_dev = true)
193  { return mfem::Write(I, Height()+1, on_dev); }
194  int *ReadWriteI(bool on_dev = true)
195  { return mfem::ReadWrite(I, Height()+1, on_dev); }
196  const int *HostReadI() const
197  { return mfem::Read(I, Height()+1, false); }
198  int *HostWriteI()
199  { return mfem::Write(I, Height()+1, false); }
201  { return mfem::ReadWrite(I, Height()+1, false); }
202 
203  // Memory access methods for the #J array.
204  Memory<int> &GetMemoryJ() { return J; }
205  const Memory<int> &GetMemoryJ() const { return J; }
206  const int *ReadJ(bool on_dev = true) const
207  { return mfem::Read(J, J.Capacity(), on_dev); }
208  int *WriteJ(bool on_dev = true)
209  { return mfem::Write(J, J.Capacity(), on_dev); }
210  int *ReadWriteJ(bool on_dev = true)
211  { return mfem::ReadWrite(J, J.Capacity(), on_dev); }
212  const int *HostReadJ() const
213  { return mfem::Read(J, J.Capacity(), false); }
214  int *HostWriteJ()
215  { return mfem::Write(J, J.Capacity(), false); }
217  { return mfem::ReadWrite(J, J.Capacity(), false); }
218 
219  // Memory access methods for the #A array.
221  const Memory<double> &GetMemoryData() const { return A; }
222  const double *ReadData(bool on_dev = true) const
223  { return mfem::Read(A, A.Capacity(), on_dev); }
224  double *WriteData(bool on_dev = true)
225  { return mfem::Write(A, A.Capacity(), on_dev); }
226  double *ReadWriteData(bool on_dev = true)
227  { return mfem::ReadWrite(A, A.Capacity(), on_dev); }
228  const double *HostReadData() const
229  { return mfem::Read(A, A.Capacity(), false); }
230  double *HostWriteData()
231  { return mfem::Write(A, A.Capacity(), false); }
233  { return mfem::ReadWrite(A, A.Capacity(), false); }
234 
235  /// Returns the number of elements in row @a i.
236  int RowSize(const int i) const;
237 
238  /// Returns the maximum number of elements among all rows.
239  int MaxRowSize() const;
240 
241  /// Return a pointer to the column indices in a row.
242  int *GetRowColumns(const int row);
243  /// Return a pointer to the column indices in a row, const version.
244  const int *GetRowColumns(const int row) const;
245 
246  /// Return a pointer to the entries in a row.
247  double *GetRowEntries(const int row);
248  /// Return a pointer to the entries in a row, const version.
249  const double *GetRowEntries(const int row) const;
250 
251  /// Change the width of a SparseMatrix.
252  /*!
253  * If width_ = -1 (DEFAULT), this routine will set the new width
254  * to the actual Width of the matrix awidth = max(J) + 1.
255  * Values 0 <= width_ < awidth are not allowed (error check in Debug Mode only)
256  *
257  * This method can be called for matrices finalized or not.
258  */
259  void SetWidth(int width_ = -1);
260 
261  /// Returns the actual Width of the matrix.
262  /*! This method can be called for matrices finalized or not. */
263  int ActualWidth() const;
264 
265  /// Sort the column indices corresponding to each row.
266  void SortColumnIndices();
267 
268  /** @brief Move the diagonal entry to the first position in each row,
269  preserving the order of the rest of the columns. */
270  void MoveDiagonalFirst();
271 
272  /// Returns reference to a_{ij}.
273  virtual double &Elem(int i, int j);
274 
275  /// Returns constant reference to a_{ij}.
276  virtual const double &Elem(int i, int j) const;
277 
278  /// Returns reference to A[i][j].
279  double &operator()(int i, int j);
280 
281  /// Returns reference to A[i][j].
282  const double &operator()(int i, int j) const;
283 
284  /// Returns the Diagonal of A
285  void GetDiag(Vector & d) const;
286 
287  /// Produces a DenseMatrix from a SparseMatrix
288  DenseMatrix *ToDenseMatrix() const;
289 
290  /// Produces a DenseMatrix from a SparseMatrix
291  void ToDenseMatrix(DenseMatrix & B) const;
292 
293  virtual MemoryClass GetMemoryClass() const
294  {
295  return Finalized() ?
297  }
298 
299  /// Matrix vector multiplication.
300  virtual void Mult(const Vector &x, Vector &y) const;
301 
302  /// y += A * x (default) or y += a * A * x
303  void AddMult(const Vector &x, Vector &y, const double a = 1.0) const;
304 
305  /// Multiply a vector with the transposed matrix. y = At * x
306  void MultTranspose(const Vector &x, Vector &y) const;
307 
308  /// y += At * x (default) or y += a * At * x
309  void AddMultTranspose(const Vector &x, Vector &y,
310  const double a = 1.0) const;
311 
312  /** @brief Build and store internally the transpose of this matrix which will
313  be used in the methods AddMultTranspose() and MultTranspose(). */
314  /** If this method has been called, the internal transpose matrix will be
315  used to perform the action of the transpose matrix in AddMultTranspose(),
316  and MultTranspose().
317 
318  Warning: any changes in this matrix will invalidate the internal
319  transpose. To rebuild the transpose, call ResetTranspose() followed by a
320  call to this method. If the internal transpose is already built, this
321  method has no effect.
322 
323  When any non-default backend is enabled, i.e. Device::IsEnabled() is
324  true, the methods AddMultTranspose(), and MultTranspose(), require the
325  internal transpose to be built. If that is not the case (i.e. the
326  internal transpose is not built), these methods will raise an error with
327  an appropriate message pointing to this method. When using the default
328  backend, calling this method is optional.
329 
330  This method can only be used when the sparse matrix is finalized. */
331  void BuildTranspose() const;
332 
333  /** Reset (destroy) the internal transpose matrix. See BuildTranspose() for
334  more details. */
335  void ResetTranspose() const;
336 
337  void PartMult(const Array<int> &rows, const Vector &x, Vector &y) const;
338  void PartAddMult(const Array<int> &rows, const Vector &x, Vector &y,
339  const double a=1.0) const;
340 
341  /// y = A * x, treating all entries as booleans (zero=false, nonzero=true).
342  /** The actual values stored in the data array, #A, are not used - this means
343  that all entries in the sparsity pattern are considered to be true by
344  this method. */
345  void BooleanMult(const Array<int> &x, Array<int> &y) const;
346 
347  /// y = At * x, treating all entries as booleans (zero=false, nonzero=true).
348  /** The actual values stored in the data array, #A, are not used - this means
349  that all entries in the sparsity pattern are considered to be true by
350  this method. */
351  void BooleanMultTranspose(const Array<int> &x, Array<int> &y) const;
352 
353  /// y = |A| * x, using entry-wise absolute values of matrix A
354  void AbsMult(const Vector &x, Vector &y) const;
355 
356  /// y = |At| * x, using entry-wise absolute values of the transpose of matrix A
357  void AbsMultTranspose(const Vector &x, Vector &y) const;
358 
359  /// Compute y^t A x
360  double InnerProduct(const Vector &x, const Vector &y) const;
361 
362  /// For all i compute \f$ x_i = \sum_j A_{ij} \f$
363  void GetRowSums(Vector &x) const;
364  /// For i = irow compute \f$ x_i = \sum_j | A_{i, j} | \f$
365  double GetRowNorml1(int irow) const;
366 
367  /// This virtual method is not supported: it always returns NULL.
368  virtual MatrixInverse *Inverse() const;
369 
370  /// Eliminates a column from the transpose matrix.
371  void EliminateRow(int row, const double sol, Vector &rhs);
372 
373  /// Eliminates a row from the matrix.
374  /*!
375  * - If @a dpolicy = #DIAG_ZERO, all the entries in the row will be set to 0.
376  * - If @a dpolicy = #DIAG_ONE (matrix must be square), the diagonal entry
377  * will be set equal to 1 and all other entries in the row to 0.
378  * - The policy #DIAG_KEEP is not supported.
379  */
380  void EliminateRow(int row, DiagonalPolicy dpolicy = DIAG_ZERO);
381 
382  /// Eliminates the column @a col from the matrix.
383  /** - If @a dpolicy = #DIAG_ZERO, all entries in the column will be set to 0.
384  - If @a dpolicy = #DIAG_ONE (matrix must be square), the diagonal entry
385  will be set equal to 1 and all other entries in the column to 0.
386  - The policy #DIAG_KEEP is not supported. */
387  void EliminateCol(int col, DiagonalPolicy dpolicy = DIAG_ZERO);
388 
389  /// Eliminate all columns i for which @a cols[i] != 0.
390  /** Elimination of a column means that all entries in the column are set to
391  zero. In addition, if the pointers @a x and @a b are not NULL, the
392  eliminated matrix entries are multiplied by the corresponding solution
393  value in @a *x and subtracted from the r.h.s. vector, @a *b. */
394  void EliminateCols(const Array<int> &cols, const Vector *x = NULL,
395  Vector *b = NULL);
396 
397  /** @brief Similar to EliminateCols + save the eliminated entries into
398  @a Ae so that (*this) + Ae is equal to the original matrix. */
399  void EliminateCols(const Array<int> &col_marker, SparseMatrix &Ae);
400 
401  /// Eliminate row @a rc and column @a rc and modify the @a rhs using @a sol.
402  /** Eliminates the column @a rc to the @a rhs, deletes the row @a rc and
403  replaces the element (rc,rc) with 1.0; assumes that element (i,rc)
404  is assembled if and only if the element (rc,i) is assembled.
405  By default, elements (rc,rc) are set to 1.0, although this behavior
406  can be adjusted by changing the @a dpolicy parameter. */
407  void EliminateRowCol(int rc, const double sol, Vector &rhs,
408  DiagonalPolicy dpolicy = DIAG_ONE);
409 
410  /** @brief Similar to
411  EliminateRowCol(int, const double, Vector &, DiagonalPolicy), but
412  multiple values for eliminated unknowns are accepted, and accordingly
413  multiple right-hand-sides are used. */
414  void EliminateRowColMultipleRHS(int rc, const Vector &sol,
415  DenseMatrix &rhs,
416  DiagonalPolicy dpolicy = DIAG_ONE);
417 
418  /// Perform elimination and set the diagonal entry to the given value
419  void EliminateRowColDiag(int rc, double value);
420 
421  /// Eliminate row @a rc and column @a rc.
422  void EliminateRowCol(int rc, DiagonalPolicy dpolicy = DIAG_ONE);
423 
424  /** @brief Similar to EliminateRowCol(int, DiagonalPolicy) + save the
425  eliminated entries into @a Ae so that (*this) + Ae is equal to the
426  original matrix */
427  void EliminateRowCol(int rc, SparseMatrix &Ae,
428  DiagonalPolicy dpolicy = DIAG_ONE);
429 
430  /// If a row contains only one diag entry of zero, set it to 1.
431  void SetDiagIdentity();
432  /// If a row contains only zeros, set its diagonal to 1.
433  virtual void EliminateZeroRows(const double threshold = 1e-12);
434 
435  /// Gauss-Seidel forward and backward iterations over a vector x.
436  void Gauss_Seidel_forw(const Vector &x, Vector &y) const;
437  void Gauss_Seidel_back(const Vector &x, Vector &y) const;
438 
439  /// Determine appropriate scaling for Jacobi iteration
440  double GetJacobiScaling() const;
441  /** One scaled Jacobi iteration for the system A x = b.
442  x1 = x0 + sc D^{-1} (b - A x0) where D is the diag of A. */
443  void Jacobi(const Vector &b, const Vector &x0, Vector &x1, double sc) const;
444 
445  void DiagScale(const Vector &b, Vector &x, double sc = 1.0) const;
446 
447  /** x1 = x0 + sc D^{-1} (b - A x0) where \f$ D_{ii} = \sum_j |A_{ij}| \f$. */
448  void Jacobi2(const Vector &b, const Vector &x0, Vector &x1,
449  double sc = 1.0) const;
450 
451  /** x1 = x0 + sc D^{-1} (b - A x0) where \f$ D_{ii} = \sum_j A_{ij} \f$. */
452  void Jacobi3(const Vector &b, const Vector &x0, Vector &x1,
453  double sc = 1.0) const;
454 
455  /** @brief Finalize the matrix initialization, switching the storage format
456  from LIL to CSR. */
457  /** This method should be called once, after the matrix has been initialized.
458  Internally, this method converts the matrix from row-wise linked list
459  (LIL) format into CSR (compressed sparse row) format. */
460  virtual void Finalize(int skip_zeros = 1) { Finalize(skip_zeros, false); }
461 
462  /// A slightly more general version of the Finalize(int) method.
463  void Finalize(int skip_zeros, bool fix_empty_rows);
464 
465  /// Returns whether or not CSR format has been finalized.
466  bool Finalized() const { return !A.Empty(); }
467  /// Returns whether or not the columns are sorted.
468  bool ColumnsAreSorted() const { return isSorted; }
469 
470  /** @brief Remove entries smaller in absolute value than a given tolerance
471  @a tol. If @a fix_empty_rows is true, a zero value is inserted in the
472  diagonal entry (for square matrices only) */
473  void Threshold(double tol, bool fix_empty_rows = false);
474 
475  /** Split the matrix into M x N blocks of sparse matrices in CSR format.
476  The 'blocks' array is M x N (i.e. M and N are determined by its
477  dimensions) and its entries are overwritten by the new blocks. */
478  void GetBlocks(Array2D<SparseMatrix *> &blocks) const;
479 
480  void GetSubMatrix(const Array<int> &rows, const Array<int> &cols,
481  DenseMatrix &subm) const;
482 
483  /** @brief Initialize the SparseMatrix for fast access to the entries of the
484  given @a row which becomes the "current row". */
485  /** Fast access to the entries of the "current row" can be performed using
486  the methods: SearchRow(const int), _Add_(const int, const double),
487  _Set_(const int, const double), and _Get_(const int). */
488  inline void SetColPtr(const int row) const;
489  /** @brief Reset the "current row" set by calling SetColPtr(). This method
490  must be called between any two calls to SetColPtr(). */
491  inline void ClearColPtr() const;
492  /// Perform a fast search for an entry in the "current row". See SetColPtr().
493  /** If the matrix is not finalized and the entry is not found in the
494  SparseMatrix, it will be added to the sparsity pattern initialized with
495  zero. If the matrix is finalized and the entry is not found, an error
496  will be generated. */
497  inline double &SearchRow(const int col);
498  /// Add a value to an entry in the "current row". See SetColPtr().
499  inline void _Add_(const int col, const double a)
500  { SearchRow(col) += a; }
501  /// Set an entry in the "current row". See SetColPtr().
502  inline void _Set_(const int col, const double a)
503  { SearchRow(col) = a; }
504  /// Read the value of an entry in the "current row". See SetColPtr().
505  inline double _Get_(const int col) const;
506 
507  inline double &SearchRow(const int row, const int col);
508  inline void _Add_(const int row, const int col, const double a)
509  { SearchRow(row, col) += a; }
510  inline void _Set_(const int row, const int col, const double a)
511  { SearchRow(row, col) = a; }
512 
513  void Set(const int i, const int j, const double a);
514  void Add(const int i, const int j, const double a);
515 
516  void SetSubMatrix(const Array<int> &rows, const Array<int> &cols,
517  const DenseMatrix &subm, int skip_zeros = 1);
518 
519  void SetSubMatrixTranspose(const Array<int> &rows, const Array<int> &cols,
520  const DenseMatrix &subm, int skip_zeros = 1);
521 
522  void AddSubMatrix(const Array<int> &rows, const Array<int> &cols,
523  const DenseMatrix &subm, int skip_zeros = 1);
524 
525  bool RowIsEmpty(const int row) const;
526 
527  /// Extract all column indices and values from a given row.
528  /** If the matrix is finalized (i.e. in CSR format), @a cols and @a srow will
529  simply be references to the specific portion of the #J and #A arrays.
530  As required by the AbstractSparseMatrix interface this method returns:
531  - 0, if @a cols and @a srow are copies of the values in the matrix, i.e.
532  when the matrix is open.
533  - 1, if @a cols and @a srow are views of the values in the matrix, i.e.
534  when the matrix is finalized.
535  @warning This method breaks the const-ness when the matrix is finalized
536  because it gives write access to the #J and #A arrays. */
537  virtual int GetRow(const int row, Array<int> &cols, Vector &srow) const;
538 
539  void SetRow(const int row, const Array<int> &cols, const Vector &srow);
540  void AddRow(const int row, const Array<int> &cols, const Vector &srow);
541 
542  void ScaleRow(const int row, const double scale);
543  /// this = diag(sl) * this;
544  void ScaleRows(const Vector & sl);
545  /// this = this * diag(sr);
546  void ScaleColumns(const Vector & sr);
547 
548  /** @brief Add the sparse matrix 'B' to '*this'. This operation will cause an
549  error if '*this' is finalized and 'B' has larger sparsity pattern. */
551 
552  /** @brief Add the sparse matrix 'B' scaled by the scalar 'a' into '*this'.
553  Only entries in the sparsity pattern of '*this' are added. */
554  void Add(const double a, const SparseMatrix &B);
555 
556  SparseMatrix &operator=(double a);
557 
558  SparseMatrix &operator*=(double a);
559 
560  /// Prints matrix to stream out.
561  void Print(std::ostream &out = mfem::out, int width_ = 4) const;
562 
563  /// Prints matrix in matlab format.
564  void PrintMatlab(std::ostream &out = mfem::out) const;
565 
566  /// Prints matrix in Matrix Market sparse format.
567  void PrintMM(std::ostream &out = mfem::out) const;
568 
569  /// Prints matrix to stream out in hypre_CSRMatrix format.
570  void PrintCSR(std::ostream &out) const;
571 
572  /// Prints a sparse matrix to stream out in CSR format.
573  void PrintCSR2(std::ostream &out) const;
574 
575  /// Print various sparse matrix statistics.
576  void PrintInfo(std::ostream &out) const;
577 
578  /// Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix
579  double IsSymmetric() const;
580 
581  /// (*this) = 1/2 ((*this) + (*this)^t)
582  void Symmetrize();
583 
584  /// Returns the number of the nonzero elements in the matrix
585  virtual int NumNonZeroElems() const;
586 
587  double MaxNorm() const;
588 
589  /// Count the number of entries with |a_ij| <= tol.
590  int CountSmallElems(double tol) const;
591 
592  /// Count the number of entries that are NOT finite, i.e. Inf or Nan.
593  int CheckFinite() const;
594 
595  /// Set the graph ownership flag (I and J arrays).
596  void SetGraphOwner(bool ownij)
597  { I.SetHostPtrOwner(ownij); J.SetHostPtrOwner(ownij); }
598 
599  /// Set the data ownership flag (A array).
600  void SetDataOwner(bool owna) { A.SetHostPtrOwner(owna); }
601 
602  /// Get the graph ownership flag (I and J arrays).
603  bool OwnsGraph() const { return I.OwnsHostPtr() && J.OwnsHostPtr(); }
604 
605  /// Get the data ownership flag (A array).
606  bool OwnsData() const { return A.OwnsHostPtr(); }
607 
608  /// Lose the ownership of the graph (I, J) and data (A) arrays.
609  void LoseData() { SetGraphOwner(false); SetDataOwner(false); }
610 
611  void Swap(SparseMatrix &other);
612 
613  /// Destroys sparse matrix.
614  virtual ~SparseMatrix()
615  {
616  Destroy();
617 #ifdef MFEM_USE_CUDA
619  {
620  cusparseDestroy(handle);
622  }
624 #endif
625  }
626 
627  Type GetType() const { return MFEM_SPARSEMAT; }
628 };
629 
630 inline std::ostream& operator<<(std::ostream& os, SparseMatrix const& mat)
631 {
632  mat.Print(os);
633  return os;
634 }
635 
636 /// Applies f() to each element of the matrix (after it is finalized).
637 void SparseMatrixFunction(SparseMatrix &S, double (*f)(double));
638 
639 
640 /// Transpose of a sparse matrix. A must be finalized.
641 SparseMatrix *Transpose(const SparseMatrix &A);
642 /// Transpose of a sparse matrix. A does not need to be a CSR matrix.
643 SparseMatrix *TransposeAbstractSparseMatrix (const AbstractSparseMatrix &A,
644  int useActualWidth);
645 
646 /// Matrix product A.B.
647 /** If @a OAB is not NULL, we assume it has the structure of A.B and store the
648  result in @a OAB. If @a OAB is NULL, we create a new SparseMatrix to store
649  the result and return a pointer to it.
650 
651  All matrices must be finalized. */
652 SparseMatrix *Mult(const SparseMatrix &A, const SparseMatrix &B,
653  SparseMatrix *OAB = NULL);
654 
655 /// C = A^T B
656 SparseMatrix *TransposeMult(const SparseMatrix &A, const SparseMatrix &B);
657 
658 /// Matrix product of sparse matrices. A and B do not need to be CSR matrices
659 SparseMatrix *MultAbstractSparseMatrix (const AbstractSparseMatrix &A,
660  const AbstractSparseMatrix &B);
661 
662 /// Matrix product A.B
663 DenseMatrix *Mult(const SparseMatrix &A, DenseMatrix &B);
664 
665 /// RAP matrix product (with R=P^T)
666 DenseMatrix *RAP(const SparseMatrix &A, DenseMatrix &P);
667 
668 /// RAP matrix product (with R=P^T)
669 DenseMatrix *RAP(DenseMatrix &A, const SparseMatrix &P);
670 
671 /** RAP matrix product (with P=R^T). ORAP is like OAB above.
672  All matrices must be finalized. */
673 SparseMatrix *RAP(const SparseMatrix &A, const SparseMatrix &R,
674  SparseMatrix *ORAP = NULL);
675 
676 /// General RAP with given R^T, A and P
677 SparseMatrix *RAP(const SparseMatrix &Rt, const SparseMatrix &A,
678  const SparseMatrix &P);
679 
680 /// Matrix multiplication A^t D A. All matrices must be finalized.
681 SparseMatrix *Mult_AtDA(const SparseMatrix &A, const Vector &D,
682  SparseMatrix *OAtDA = NULL);
683 
684 
685 /// Matrix addition result = A + B.
686 SparseMatrix * Add(const SparseMatrix & A, const SparseMatrix & B);
687 /// Matrix addition result = a*A + b*B
688 SparseMatrix * Add(double a, const SparseMatrix & A, double b,
689  const SparseMatrix & B);
690 /// Matrix addition result = sum_i A_i
691 SparseMatrix * Add(Array<SparseMatrix *> & Ai);
692 
693 /// B += alpha * A
694 void Add(const SparseMatrix &A, double alpha, DenseMatrix &B);
695 
696 /// Produces a block matrix with blocks A_{ij}*B
697 DenseMatrix *OuterProduct(const DenseMatrix &A, const DenseMatrix &B);
698 
699 /// Produces a block matrix with blocks A_{ij}*B
700 SparseMatrix *OuterProduct(const DenseMatrix &A, const SparseMatrix &B);
701 
702 /// Produces a block matrix with blocks A_{ij}*B
703 SparseMatrix *OuterProduct(const SparseMatrix &A, const DenseMatrix &B);
704 
705 /// Produces a block matrix with blocks A_{ij}*B
706 SparseMatrix *OuterProduct(const SparseMatrix &A, const SparseMatrix &B);
707 
708 
709 // Inline methods
710 
711 inline void SparseMatrix::SetColPtr(const int row) const
712 {
713  if (Rows)
714  {
715  if (ColPtrNode == NULL)
716  {
717  ColPtrNode = new RowNode *[width];
718  for (int i = 0; i < width; i++)
719  {
720  ColPtrNode[i] = NULL;
721  }
722  }
723  for (RowNode *node_p = Rows[row]; node_p != NULL; node_p = node_p->Prev)
724  {
725  ColPtrNode[node_p->Column] = node_p;
726  }
727  }
728  else
729  {
730  if (ColPtrJ == NULL)
731  {
732  ColPtrJ = new int[width];
733  for (int i = 0; i < width; i++)
734  {
735  ColPtrJ[i] = -1;
736  }
737  }
738  for (int j = I[row], end = I[row+1]; j < end; j++)
739  {
740  ColPtrJ[J[j]] = j;
741  }
742  }
743  current_row = row;
744 }
745 
746 inline void SparseMatrix::ClearColPtr() const
747 {
748  if (Rows)
749  {
750  for (RowNode *node_p = Rows[current_row]; node_p != NULL;
751  node_p = node_p->Prev)
752  {
753  ColPtrNode[node_p->Column] = NULL;
754  }
755  }
756  else
757  {
758  for (int j = I[current_row], end = I[current_row+1]; j < end; j++)
759  {
760  ColPtrJ[J[j]] = -1;
761  }
762  }
763 }
764 
765 inline double &SparseMatrix::SearchRow(const int col)
766 {
767  if (Rows)
768  {
769  RowNode *node_p = ColPtrNode[col];
770  if (node_p == NULL)
771  {
772 #ifdef MFEM_USE_MEMALLOC
773  node_p = NodesMem->Alloc();
774 #else
775  node_p = new RowNode;
776 #endif
777  node_p->Prev = Rows[current_row];
778  node_p->Column = col;
779  node_p->Value = 0.0;
780  Rows[current_row] = ColPtrNode[col] = node_p;
781  }
782  return node_p->Value;
783  }
784  else
785  {
786  const int j = ColPtrJ[col];
787  MFEM_VERIFY(j != -1, "Entry for column " << col << " is not allocated.");
788  return A[j];
789  }
790 }
791 
792 inline double SparseMatrix::_Get_(const int col) const
793 {
794  if (Rows)
795  {
796  RowNode *node_p = ColPtrNode[col];
797  return (node_p == NULL) ? 0.0 : node_p->Value;
798  }
799  else
800  {
801  const int j = ColPtrJ[col];
802  return (j == -1) ? 0.0 : A[j];
803  }
804 }
805 
806 inline double &SparseMatrix::SearchRow(const int row, const int col)
807 {
808  if (Rows)
809  {
810  RowNode *node_p;
811 
812  for (node_p = Rows[row]; 1; node_p = node_p->Prev)
813  {
814  if (node_p == NULL)
815  {
816 #ifdef MFEM_USE_MEMALLOC
817  node_p = NodesMem->Alloc();
818 #else
819  node_p = new RowNode;
820 #endif
821  node_p->Prev = Rows[row];
822  node_p->Column = col;
823  node_p->Value = 0.0;
824  Rows[row] = node_p;
825  break;
826  }
827  else if (node_p->Column == col)
828  {
829  break;
830  }
831  }
832  return node_p->Value;
833  }
834  else
835  {
836  int *Ip = I+row, *Jp = J;
837  for (int k = Ip[0], end = Ip[1]; k < end; k++)
838  {
839  if (Jp[k] == col)
840  {
841  return A[k];
842  }
843  }
844  MFEM_ABORT("Could not find entry for row = " << row << ", col = " << col);
845  }
846  return A[0];
847 }
848 
849 /// Specialization of the template function Swap<> for class SparseMatrix
850 template<> inline void Swap<SparseMatrix>(SparseMatrix &a, SparseMatrix &b)
851 {
852  a.Swap(b);
853 }
854 
855 } // namespace mfem
856 
857 #endif
RowNode ** ColPtrNode
Definition: sparsemat.hpp:73
Memory< int > I
Array with size (height+1) containing the row offsets.
Definition: sparsemat.hpp:58
Elem * Alloc()
Definition: mem_alloc.hpp:166
void _Add_(const int row, const int col, const double a)
Definition: sparsemat.hpp:508
double GetJacobiScaling() const
Determine appropriate scaling for Jacobi iteration.
Definition: sparsemat.cpp:2277
int CheckFinite() const
Count the number of entries that are NOT finite, i.e. Inf or Nan.
Definition: sparsemat.cpp:1434
cusparseStatus_t status
Definition: sparsemat.hpp:95
SparseMatrix * TransposeAbstractSparseMatrix(const AbstractSparseMatrix &A, int useActualWidth)
Transpose of a sparse matrix. A does not need to be a CSR matrix.
Definition: sparsemat.cpp:3262
int RowSize(const int i) const
Returns the number of elements in row i.
Definition: sparsemat.cpp:296
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:1359
void Add(const int i, const int j, const double a)
Definition: sparsemat.cpp:2493
void * CuMemFree(void *dptr)
Frees device memory and returns destination ptr.
Definition: cuda.cpp:64
void _Add_(const int col, const double a)
Add a value to an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:499
void BuildTranspose() const
Build and store internally the transpose of this matrix which will be used in the methods AddMultTran...
Definition: sparsemat.cpp:773
void Clear()
Clear the contents of the SparseMatrix.
Definition: sparsemat.hpp:167
void EliminateCol(int col, DiagonalPolicy dpolicy=DIAG_ZERO)
Eliminates the column col from the matrix.
Definition: sparsemat.cpp:1511
SparseMatrix * At
Transpose of A. Owned. Used to perform MultTranspose() on devices.
Definition: sparsemat.hpp:76
Type GetType() const
Definition: sparsemat.hpp:627
static MemoryClass GetHostMemoryClass()
Get the current Host MemoryClass. This is the MemoryClass used by most MFEM host Memory objects...
Definition: device.hpp:256
bool ColumnsAreSorted() const
Returns whether or not the columns are sorted.
Definition: sparsemat.hpp:468
void _Set_(const int row, const int col, const double a)
Definition: sparsemat.hpp:510
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1640
static size_t bufferSize
Definition: sparsemat.hpp:98
void DiagScale(const Vector &b, Vector &x, double sc=1.0) const
Definition: sparsemat.cpp:2341
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
Definition: sparsemat.hpp:460
cusparseMatDescr_t descr
Definition: sparsemat.hpp:97
void MakeRef(const SparseMatrix &master)
Clear the contents of the SparseMatrix and make it a reference to master.
Definition: sparsemat.cpp:257
double & SearchRow(const int col)
Perform a fast search for an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:765
DenseMatrix * ToDenseMatrix() const
Produces a DenseMatrix from a SparseMatrix.
Definition: sparsemat.cpp:564
bool Empty() const
Check if the SparseMatrix is empty.
Definition: sparsemat.hpp:170
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Definition: sparsemat.cpp:343
void BooleanMult(const Array< int > &x, Array< int > &y) const
y = A * x, treating all entries as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:831
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
const double * ReadData(bool on_dev=true) const
Definition: sparsemat.hpp:222
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:178
Abstract data type for sparse matrices.
Definition: matrix.hpp:73
SparseMatrix & operator+=(const SparseMatrix &B)
Add the sparse matrix &#39;B&#39; to &#39;*this&#39;. This operation will cause an error if &#39;*this&#39; is finalized and ...
Definition: sparsemat.cpp:2835
T * Write(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for write access to mem with the mfem::Device&#39;s DeviceMemoryClass, if on_dev = true...
Definition: device.hpp:324
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
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:601
double MaxNorm() const
Definition: sparsemat.cpp:1381
virtual ~SparseMatrix()
Destroys sparse matrix.
Definition: sparsemat.hpp:614
int * GetI()
Return the array I.
Definition: sparsemat.hpp:173
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:728
Abstract data type for matrix inverse.
Definition: matrix.hpp:62
static void * dBuffer
Definition: sparsemat.hpp:99
void EliminateRowColMultipleRHS(int rc, const Vector &sol, DenseMatrix &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Similar to EliminateRowCol(int, const double, Vector &amp;, DiagonalPolicy), but multiple values for elim...
Definition: sparsemat.cpp:1720
void PrintInfo(std::ostream &out) const
Print various sparse matrix statistics.
Definition: sparsemat.cpp:3077
void UseCuSparse(bool _useCuSparse=true)
Definition: sparsemat.hpp:152
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
Definition: sparsemat.cpp:462
void LoseData()
Lose the ownership of the graph (I, J) and data (A) arrays.
Definition: sparsemat.hpp:609
void AddRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:2716
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
Definition: sparsemat.cpp:357
void SortColumnIndices()
Sort the column indices corresponding to each row.
Definition: sparsemat.cpp:409
Memory< double > & GetMemoryData()
Definition: sparsemat.hpp:220
void PrintMM(std::ostream &out=mfem::out) const
Prints matrix in Matrix Market sparse format.
Definition: sparsemat.cpp:3006
int Capacity() const
Return the size of the allocated memory.
Memory< int > & GetMemoryI()
Definition: sparsemat.hpp:188
cusparseDnVecDescr_t vecY_descr
Definition: sparsemat.hpp:105
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const
Definition: sparsemat.cpp:2580
double * HostReadWriteData()
Definition: sparsemat.hpp:232
void _Set_(const int col, const double a)
Set an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:502
void BooleanMultTranspose(const Array< int > &x, Array< int > &y) const
y = At * x, treating all entries as booleans (zero=false, nonzero=true).
Definition: sparsemat.cpp:861
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:3630
void PartAddMult(const Array< int > &rows, const Vector &x, Vector &y, const double a=1.0) const
Definition: sparsemat.cpp:813
bool RowIsEmpty(const int row) const
Definition: sparsemat.cpp:2608
void GetBlocks(Array2D< SparseMatrix * > &blocks) const
Definition: sparsemat.cpp:1208
const int * ReadJ(bool on_dev=true) const
Definition: sparsemat.hpp:206
void ScaleRow(const int row, const double scale)
Definition: sparsemat.cpp:2748
RowNodeAlloc * NodesMem
Definition: sparsemat.hpp:80
double * GetData()
Return the element data, i.e. the array A.
Definition: sparsemat.hpp:183
Memory< int > & GetMemoryJ()
Definition: sparsemat.hpp:204
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
Definition: sparsemat.cpp:1340
void MoveDiagonalFirst()
Move the diagonal entry to the first position in each row, preserving the order of the rest of the co...
Definition: sparsemat.cpp:438
void SetGraphOwner(bool ownij)
Set the graph ownership flag (I and J arrays).
Definition: sparsemat.hpp:596
ID for class SparseMatrix.
Definition: operator.hpp:240
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:1928
bool OwnsData() const
Get the data ownership flag (A array).
Definition: sparsemat.hpp:606
void ClearColPtr() const
Reset the &quot;current row&quot; set by calling SetColPtr(). This method must be called between any two calls ...
Definition: sparsemat.hpp:746
void ScaleColumns(const Vector &sr)
this = this * diag(sr);
Definition: sparsemat.cpp:2807
void PrintCSR2(std::ostream &out) const
Prints a sparse matrix to stream out in CSR format.
Definition: sparsemat.cpp:3052
virtual int GetRow(const int row, Array< int > &cols, Vector &srow) const
Extract all column indices and values from a given row.
Definition: sparsemat.cpp:2629
bool isSorted
Are the columns sorted already.
Definition: sparsemat.hpp:84
Memory< double > A
Array with size I[height], containing the actual entries of the sparse matrix, as indexed by the I ar...
Definition: sparsemat.hpp:64
SparseMatrix()
Create an empty SparseMatrix.
Definition: sparsemat.hpp:110
Data type sparse matrix.
Definition: sparsemat.hpp:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
static MemoryClass GetDeviceMemoryClass()
Get the current Device MemoryClass. This is the MemoryClass used by most MFEM device kernels to acces...
Definition: device.hpp:269
bool OwnsGraph() const
Get the graph ownership flag (I and J arrays).
Definition: sparsemat.hpp:603
void SetWidth(int width_=-1)
Change the width of a SparseMatrix.
Definition: sparsemat.cpp:371
double b
Definition: lissajous.cpp:42
const double * HostReadData() const
Definition: sparsemat.hpp:228
int * ReadWriteI(bool on_dev=true)
Definition: sparsemat.hpp:194
double IsSymmetric() const
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix.
Definition: sparsemat.cpp:1298
const int * HostReadJ() const
Definition: sparsemat.hpp:212
static int SparseMatrixCount
Definition: sparsemat.hpp:102
void SetDiagIdentity()
If a row contains only one diag entry of zero, set it to 1.
Definition: sparsemat.cpp:2080
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
Definition: table.cpp:414
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;
Definition: sparsemat.cpp:2776
static cusparseHandle_t handle
Definition: sparsemat.hpp:96
void Swap< SparseMatrix >(SparseMatrix &a, SparseMatrix &b)
Specialization of the template function Swap&lt;&gt; for class SparseMatrix.
Definition: sparsemat.hpp:850
DenseMatrix * OuterProduct(const DenseMatrix &A, const DenseMatrix &B)
Produces a block matrix with blocks A_{ij}*B.
Definition: sparsemat.cpp:3774
int * WriteI(bool on_dev=true)
Definition: sparsemat.hpp:192
Set the diagonal value to one.
Definition: operator.hpp:48
void PrintMatlab(std::ostream &out=mfem::out) const
Prints matrix in matlab format.
Definition: sparsemat.cpp:2986
const T * Read(const Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read access to mem with the mfem::Device&#39;s DeviceMemoryClass, if on_dev = true...
Definition: device.hpp:307
Biwise-OR of all CUDA backends.
Definition: device.hpp:85
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2436
Dynamic 2D array using row-major layout.
Definition: array.hpp:333
void Jacobi(const Vector &b, const Vector &x0, Vector &x1, double sc) const
Definition: sparsemat.cpp:2310
void SetSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2512
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
Definition: sparsemat.hpp:293
bool Finalized() const
Returns whether or not CSR format has been finalized.
Definition: sparsemat.hpp:466
const double * GetData() const
Return the element data, i.e. the array A, const version.
Definition: sparsemat.hpp:185
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:237
void SparseMatrixFunction(SparseMatrix &S, double(*f)(double))
Applies f() to each element of the matrix (after it is finalized).
Definition: sparsemat.cpp:3195
virtual MatrixInverse * Inverse() const
This virtual method is not supported: it always returns NULL.
Definition: sparsemat.cpp:1458
const Memory< double > & GetMemoryData() const
Definition: sparsemat.hpp:221
void EliminateRowColDiag(int rc, double value)
Perform elimination and set the diagonal entry to the given value.
Definition: sparsemat.cpp:1923
const Memory< int > & GetMemoryI() const
Definition: sparsemat.hpp:189
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:594
void SetColPtr(const int row) const
Initialize the SparseMatrix for fast access to the entries of the given row which becomes the &quot;curren...
Definition: sparsemat.hpp:711
const int * GetJ() const
Return the array J, const version.
Definition: sparsemat.hpp:180
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)
Definition: sparsemat.cpp:2545
void Set(const int i, const int j, const double a)
Definition: sparsemat.cpp:2474
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.
Definition: sparsemat.cpp:1624
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...
Definition: sparsemat.cpp:1067
double * HostWriteData()
Definition: sparsemat.hpp:230
int Size() const
For backward compatibility, define Size() to be synonym of Height().
Definition: sparsemat.hpp:164
void Gauss_Seidel_back(const Vector &x, Vector &y) const
Definition: sparsemat.cpp:2194
const int * ReadI(bool on_dev=true) const
Definition: sparsemat.hpp:190
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
Definition: device.hpp:246
void Gauss_Seidel_forw(const Vector &x, Vector &y) const
Gauss-Seidel forward and backward iterations over a vector x.
Definition: sparsemat.cpp:2110
int * HostReadWriteJ()
Definition: sparsemat.hpp:216
RowNode ** Rows
Array of linked lists, one for every row. This array represents the linked list (LIL) storage format...
Definition: sparsemat.hpp:69
const Memory< int > & GetMemoryJ() const
Definition: sparsemat.hpp:205
int * WriteJ(bool on_dev=true)
Definition: sparsemat.hpp:208
void EliminateRow(int row, const double sol, Vector &rhs)
Eliminates a column from the transpose matrix.
Definition: sparsemat.cpp:1463
double a
Definition: lissajous.cpp:41
int ActualWidth() const
Returns the actual Width of the matrix.
Definition: sparsemat.cpp:3169
void GetDiag(Vector &d) const
Returns the Diagonal of A.
Definition: sparsemat.cpp:530
MemAlloc< RowNode, 1024 > RowNodeAlloc
Definition: sparsemat.hpp:79
T * ReadWrite(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read+write access to mem with the mfem::Device&#39;s DeviceMemoryClass, if on_dev = true, or the mfem::Device&#39;s HostMemoryClass, otherwise.
Definition: device.hpp:341
double GetRowNorml1(int irow) const
For i = irow compute .
Definition: sparsemat.cpp:1043
void ResetTranspose() const
Definition: sparsemat.cpp:781
SparseMatrix * TransposeMult(const SparseMatrix &A, const SparseMatrix &B)
C = A^T B.
Definition: sparsemat.cpp:3467
virtual void EliminateZeroRows(const double threshold=1e-12)
If a row contains only zeros, set its diagonal to 1.
Definition: sparsemat.cpp:2091
bool Empty() const
Return true if the Memory object is empty, see Reset().
void AbsMultTranspose(const Vector &x, Vector &y) const
y = |At| * x, using entry-wise absolute values of the transpose of matrix A
Definition: sparsemat.cpp:933
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Definition: operator.hpp:45
double & operator()(int i, int j)
Returns reference to A[i][j].
Definition: sparsemat.cpp:472
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:721
void EliminateCols(const Array< int > &cols, const Vector *x=NULL, Vector *b=NULL)
Eliminate all columns i for which cols[i] != 0.
Definition: sparsemat.cpp:1551
SparseMatrix & operator=(const SparseMatrix &rhs)
Assignment operator: deep copy.
Definition: sparsemat.cpp:247
SparseMatrix & operator*=(double a)
Definition: sparsemat.cpp:2914
void SetDataOwner(bool owna)
Set the data ownership flag (A array).
Definition: sparsemat.hpp:600
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
Definition: sparsemat.hpp:630
void SetRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:2669
cusparseSpMatDescr_t matA_descr
Definition: sparsemat.hpp:103
const double alpha
Definition: ex15.cpp:336
double * ReadWriteData(bool on_dev=true)
Definition: sparsemat.hpp:226
int * ReadWriteJ(bool on_dev=true)
Definition: sparsemat.hpp:210
Vector data type.
Definition: vector.hpp:51
void Print(std::ostream &out=mfem::out, int width_=4) const
Prints matrix to stream out.
Definition: sparsemat.cpp:2938
const int * HostReadI() const
Definition: sparsemat.hpp:196
const int * GetI() const
Return the array I, const version.
Definition: sparsemat.hpp:175
void Jacobi2(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
Definition: sparsemat.cpp:2388
cusparseDnVecDescr_t vecX_descr
Definition: sparsemat.hpp:104
int MaxRowSize() const
Returns the maximum number of elements among all rows.
Definition: sparsemat.cpp:319
int * HostReadWriteI()
Definition: sparsemat.hpp:200
void Jacobi3(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
Definition: sparsemat.cpp:2412
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:3475
void AbsMult(const Vector &x, Vector &y) const
y = |A| * x, using entry-wise absolute values of matrix A
Definition: sparsemat.cpp:884
double _Get_(const int col) const
Read the value of an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:792
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
Memory< int > J
Array with size I[height], containing the column indices for all matrix entries, as indexed by the I ...
Definition: sparsemat.hpp:61
void Swap(SparseMatrix &other)
Definition: sparsemat.cpp:3882
void PrintCSR(std::ostream &out) const
Prints matrix to stream out in hypre_CSRMatrix format.
Definition: sparsemat.cpp:3028
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:979
MemoryClass
Memory classes identify sets of memory types.
Definition: mem_manager.hpp:58
Set the diagonal value to zero.
Definition: operator.hpp:47
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
int CountSmallElems(double tol) const
Count the number of entries with |a_ij| &lt;= tol.
Definition: sparsemat.cpp:1406
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const
Definition: sparsemat.cpp:787
void GetRowSums(Vector &x) const
For all i compute .
Definition: sparsemat.cpp:1020
double f(const Vector &p)
double * WriteData(bool on_dev=true)
Definition: sparsemat.hpp:224