MFEM  v4.5.1
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-2022, 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/backends.hpp"
18 #include "../general/mem_alloc.hpp"
19 #include "../general/mem_manager.hpp"
20 #include "../general/device.hpp"
21 #include "../general/table.hpp"
22 #include "../general/globals.hpp"
23 #include "densemat.hpp"
24 
25 #if defined(MFEM_USE_HIP)
26 #if (HIP_VERSION_MAJOR * 100 + HIP_VERSION_MINOR) < 502
27 #include <hipsparse.h>
28 #else
29 #include <hipsparse/hipsparse.h>
30 #endif
31 #endif
32 
33 
34 namespace mfem
35 {
36 
37 class
38 #if defined(__alignas_is_defined)
39  alignas(double)
40 #endif
41  RowNode
42 {
43 public:
44  double Value;
45  RowNode *Prev;
46  int Column;
47 };
48 
49 /// Data type sparse matrix
51 {
52 protected:
53  /// @name Arrays used by the CSR storage format.
54  /** */
55  ///@{
56  /// @brief %Array with size (#height+1) containing the row offsets.
57  /** The data for row r, 0 <= r < height, is at offsets j, I[r] <= j < I[r+1].
58  The offsets, j, are indices in the #J and #A arrays. The first entry in
59  this array is always zero, I[0] = 0, and the last entry, I[height], gives
60  the total number of entries stored (at a minimum, all nonzeros must be
61  represented) in the sparse matrix. */
63  /** @brief %Array with size #I[#height], containing the column indices for
64  all matrix entries, as indexed by the #I array. */
66  /** @brief %Array with size #I[#height], containing the actual entries of the
67  sparse matrix, as indexed by the #I array. */
69  ///@}
70 
71  /** @brief %Array of linked lists, one for every row. This array represents
72  the linked list (LIL) storage format. */
73  RowNode **Rows;
74 
75  mutable int current_row;
76  mutable int* ColPtrJ;
77  mutable RowNode ** ColPtrNode;
78 
79  /// Transpose of A. Owned. Used to perform MultTranspose() on devices.
80  mutable SparseMatrix *At;
81 
82 #ifdef MFEM_USE_MEMALLOC
85 #endif
86 
87  /// Are the columns sorted already.
88  bool isSorted;
89 
90  void Destroy(); // Delete all owned data
91  void SetEmpty(); // Init all entries with empty values
92 
93  bool useGPUSparse = true; // Use cuSPARSE or hipSPARSE if available
94 
95  // Initialize cuSPARSE/hipSPARSE
96  void InitGPUSparse();
97 
98 #ifdef MFEM_USE_CUDA_OR_HIP
99  // common for hipSPARSE and cuSPARSE
100  static int SparseMatrixCount;
101  static size_t bufferSize;
102  static void *dBuffer;
103  mutable bool initBuffers = false;
104 
105 #if defined(MFEM_USE_CUDA)
106  cusparseStatus_t status;
107  static cusparseHandle_t handle;
108  cusparseMatDescr_t descr = 0;
109 
110 #if CUDA_VERSION >= 10010
111  mutable cusparseSpMatDescr_t matA_descr;
112  mutable cusparseDnVecDescr_t vecX_descr;
113  mutable cusparseDnVecDescr_t vecY_descr;
114 #else // CUDA_VERSION >= 10010
115  mutable cusparseMatDescr_t matA_descr;
116 #endif // CUDA_VERSION >= 10010
117 
118 #else // defined(MFEM_USE_CUDA)
119  hipsparseStatus_t status;
120  static hipsparseHandle_t handle;
121  hipsparseMatDescr_t descr = 0;
122 
123  mutable hipsparseSpMatDescr_t matA_descr;
124  mutable hipsparseDnVecDescr_t vecX_descr;
125  mutable hipsparseDnVecDescr_t vecY_descr;
126 #endif // defined(MFEM_USE_CUDA)
127 #endif // MFEM_USE_CUDA_OR_HIP
128 
129 public:
130  /// Create an empty SparseMatrix.
132  {
133  SetEmpty();
134 
135  InitGPUSparse();
136  }
137 
138  /** @brief Create a sparse matrix with flexible sparsity structure using a
139  row-wise linked list (LIL) format. */
140  /** New entries are added as needed by methods like AddSubMatrix(),
141  SetSubMatrix(), etc. Calling Finalize() will convert the SparseMatrix to
142  the more compact compressed sparse row (CSR) format. */
143  explicit SparseMatrix(int nrows, int ncols = -1);
144 
145  /** @brief Create a sparse matrix in CSR format. Ownership of @a i, @a j, and
146  @a data is transferred to the SparseMatrix. */
147  SparseMatrix(int *i, int *j, double *data, int m, int n);
148 
149  /** @brief Create a sparse matrix in CSR format. Ownership of @a i, @a j, and
150  @a data is optionally transferred to the SparseMatrix. */
151  /** If the parameter @a data is NULL, then the internal #A array is allocated
152  by this constructor (initializing it with zeros and taking ownership,
153  regardless of the parameter @a owna). */
154  SparseMatrix(int *i, int *j, double *data, int m, int n, bool ownij,
155  bool owna, bool issorted);
156 
157  /** @brief Create a sparse matrix in CSR format where each row has space
158  allocated for exactly @a rowsize entries. */
159  /** SetRow() can then be called or the #I, #J, #A arrays can be used
160  directly. */
161  SparseMatrix(int nrows, int ncols, int rowsize);
162 
163  /// Copy constructor (deep copy).
164  /** If @a mat is finalized and @a copy_graph is false, the #I and #J arrays
165  will use a shallow copy (copy the pointers only) without transferring
166  ownership.
167  If @a mt is MemoryType::PRESERVE the memory type of the resulting
168  SparseMatrix's #I, #J, and #A arrays will be the same as @a mat,
169  otherwise the type will be @a mt for those arrays that are deep
170  copied. */
171  SparseMatrix(const SparseMatrix &mat, bool copy_graph = true,
173 
174  /// Create a SparseMatrix with diagonal @a v, i.e. A = Diag(v)
175  SparseMatrix(const Vector & v);
176 
177  /// @brief Sets the height and width of the matrix.
178  /** @warning This does not modify in any way the underlying CSR or LIL
179  representation of the matrix.
180 
181  This function should generally be called when manually constructing the
182  CSR #I, #J, and #A arrays in conjunction with the
183  SparseMatrix::SparseMatrix() constructor. */
184  void OverrideSize(int height_, int width_);
185 
186  /** @brief Runtime option to use cuSPARSE or hipSPARSE. Only valid when using
187  a CUDA or HIP backend.
188 
189  @note This option is enabled by default, so typically one would use this
190  method to disable the use of cuSPARSE/hipSPARSE. */
191  void UseGPUSparse(bool useGPUSparse_ = true) { useGPUSparse = useGPUSparse_;}
192  /// Deprecated equivalent of UseGPUSparse().
193  MFEM_DEPRECATED
194  void UseCuSparse(bool useCuSparse_ = true) { UseGPUSparse(useCuSparse_); }
195 
196  /// Assignment operator: deep copy
197  SparseMatrix& operator=(const SparseMatrix &rhs);
198 
199  /** @brief Clear the contents of the SparseMatrix and make it a reference to
200  @a master */
201  /** After this call, the matrix will point to the same data as @a master but
202  it will not own its data. The @a master must be finalized. */
203  void MakeRef(const SparseMatrix &master);
204 
205  /// For backward compatibility, define Size() to be synonym of Height().
206  int Size() const { return Height(); }
207 
208  /// Clear the contents of the SparseMatrix.
209  void Clear() { Destroy(); SetEmpty(); }
210 
211  /** @brief Clear the cuSPARSE/hipSPARSE descriptors.
212  This must be called after releasing the device memory of A. */
213  void ClearGPUSparse();
214  /// Deprecated equivalent of ClearGPUSparse().
215  MFEM_DEPRECATED
217 
218  /// Check if the SparseMatrix is empty.
219  bool Empty() const { return (A == NULL) && (Rows == NULL); }
220 
221  /// Return the array #I.
222  inline int *GetI() { return I; }
223  /// Return the array #I, const version.
224  inline const int *GetI() const { return I; }
225 
226  /// Return the array #J.
227  inline int *GetJ() { return J; }
228  /// Return the array #J, const version.
229  inline const int *GetJ() const { return J; }
230 
231  /// Return the element data, i.e. the array #A.
232  inline double *GetData() { return A; }
233  /// Return the element data, i.e. the array #A, const version.
234  inline const double *GetData() const { return A; }
235 
236  // Memory access methods for the #I array.
237  Memory<int> &GetMemoryI() { return I; }
238  const Memory<int> &GetMemoryI() const { return I; }
239  const int *ReadI(bool on_dev = true) const
240  { return mfem::Read(I, Height()+1, on_dev); }
241  int *WriteI(bool on_dev = true)
242  { return mfem::Write(I, Height()+1, on_dev); }
243  int *ReadWriteI(bool on_dev = true)
244  { return mfem::ReadWrite(I, Height()+1, on_dev); }
245  const int *HostReadI() const
246  { return mfem::Read(I, Height()+1, false); }
247  int *HostWriteI()
248  { return mfem::Write(I, Height()+1, false); }
250  { return mfem::ReadWrite(I, Height()+1, false); }
251 
252  // Memory access methods for the #J array.
253  Memory<int> &GetMemoryJ() { return J; }
254  const Memory<int> &GetMemoryJ() const { return J; }
255  const int *ReadJ(bool on_dev = true) const
256  { return mfem::Read(J, J.Capacity(), on_dev); }
257  int *WriteJ(bool on_dev = true)
258  { return mfem::Write(J, J.Capacity(), on_dev); }
259  int *ReadWriteJ(bool on_dev = true)
260  { return mfem::ReadWrite(J, J.Capacity(), on_dev); }
261  const int *HostReadJ() const
262  { return mfem::Read(J, J.Capacity(), false); }
263  int *HostWriteJ()
264  { return mfem::Write(J, J.Capacity(), false); }
266  { return mfem::ReadWrite(J, J.Capacity(), false); }
267 
268  // Memory access methods for the #A array.
270  const Memory<double> &GetMemoryData() const { return A; }
271  const double *ReadData(bool on_dev = true) const
272  { return mfem::Read(A, A.Capacity(), on_dev); }
273  double *WriteData(bool on_dev = true)
274  { return mfem::Write(A, A.Capacity(), on_dev); }
275  double *ReadWriteData(bool on_dev = true)
276  { return mfem::ReadWrite(A, A.Capacity(), on_dev); }
277  const double *HostReadData() const
278  { return mfem::Read(A, A.Capacity(), false); }
279  double *HostWriteData()
280  { return mfem::Write(A, A.Capacity(), false); }
282  { return mfem::ReadWrite(A, A.Capacity(), false); }
283 
284  /// Returns the number of elements in row @a i.
285  int RowSize(const int i) const;
286 
287  /// Returns the maximum number of elements among all rows.
288  int MaxRowSize() const;
289 
290  /// Return a pointer to the column indices in a row.
291  int *GetRowColumns(const int row);
292  /// Return a pointer to the column indices in a row, const version.
293  const int *GetRowColumns(const int row) const;
294 
295  /// Return a pointer to the entries in a row.
296  double *GetRowEntries(const int row);
297  /// Return a pointer to the entries in a row, const version.
298  const double *GetRowEntries(const int row) const;
299 
300  /// Change the width of a SparseMatrix.
301  /*!
302  * If width_ = -1 (DEFAULT), this routine will set the new width
303  * to the actual Width of the matrix awidth = max(J) + 1.
304  * Values 0 <= width_ < awidth are not allowed (error check in Debug Mode only)
305  *
306  * This method can be called for matrices finalized or not.
307  */
308  void SetWidth(int width_ = -1);
309 
310  /// Returns the actual Width of the matrix.
311  /*! This method can be called for matrices finalized or not. */
312  int ActualWidth() const;
313 
314  /// Sort the column indices corresponding to each row.
315  void SortColumnIndices();
316 
317  /** @brief Move the diagonal entry to the first position in each row,
318  preserving the order of the rest of the columns. */
319  void MoveDiagonalFirst();
320 
321  /// Returns reference to a_{ij}.
322  virtual double &Elem(int i, int j);
323 
324  /// Returns constant reference to a_{ij}.
325  virtual const double &Elem(int i, int j) const;
326 
327  /// Returns reference to A[i][j].
328  double &operator()(int i, int j);
329 
330  /// Returns reference to A[i][j].
331  const double &operator()(int i, int j) const;
332 
333  /// Returns the Diagonal of A
334  void GetDiag(Vector & d) const;
335 
336  /// Produces a DenseMatrix from a SparseMatrix
337  DenseMatrix *ToDenseMatrix() const;
338 
339  /// Produces a DenseMatrix from a SparseMatrix
340  void ToDenseMatrix(DenseMatrix & B) const;
341 
342  virtual MemoryClass GetMemoryClass() const
343  {
344  return Finalized() ?
346  }
347 
348  /// Matrix vector multiplication.
349  virtual void Mult(const Vector &x, Vector &y) const;
350 
351  /// y += A * x (default) or y += a * A * x
352  void AddMult(const Vector &x, Vector &y, const double a = 1.0) const;
353 
354  /// Multiply a vector with the transposed matrix. y = At * x
355  /** If the matrix is modified, call ResetTranspose() and optionally
356  EnsureMultTranspose() to make sure this method uses the correct updated
357  transpose. */
358  void MultTranspose(const Vector &x, Vector &y) const;
359 
360  /// y += At * x (default) or y += a * At * x
361  /** If the matrix is modified, call ResetTranspose() and optionally
362  EnsureMultTranspose() to make sure this method uses the correct updated
363  transpose. */
364  void AddMultTranspose(const Vector &x, Vector &y,
365  const double a = 1.0) const;
366 
367  /** @brief Build and store internally the transpose of this matrix which will
368  be used in the methods AddMultTranspose(), MultTranspose(), and
369  AbsMultTranspose(). */
370  /** If this method has been called, the internal transpose matrix will be
371  used to perform the action of the transpose matrix in AddMultTranspose(),
372  MultTranspose(), and AbsMultTranspose().
373 
374  Warning: any changes in this matrix will invalidate the internal
375  transpose. To rebuild the transpose, call ResetTranspose() followed by
376  (optionally) a call to this method. If the internal transpose is already
377  built, this method has no effect.
378 
379  When any non-serial-CPU backend is enabled, i.e. the call
380  Device::Allows(~ Backend::CPU_MASK) returns true, the above methods
381  require the internal transpose to be built. If that is not the case (i.e.
382  the internal transpose is not built), these methods will automatically
383  call EnsureMultTranspose(). When using any backend from
384  Backend::CPU_MASK, calling this method is optional.
385 
386  This method can only be used when the sparse matrix is finalized.
387 
388  @sa EnsureMultTranspose(), ResetTranspose(). */
389  void BuildTranspose() const;
390 
391  /** Reset (destroy) the internal transpose matrix. See BuildTranspose() for
392  more details. */
393  void ResetTranspose() const;
394 
395  /** @brief Ensures that the matrix is capable of performing MultTranspose(),
396  AddMultTranspose(), and AbsMultTranspose(). */
397  /** For non-serial-CPU backends (e.g. GPU, OpenMP), multiplying by the
398  transpose requires that the internal transpose matrix be already built.
399  When such a backend is enabled, this function will build the internal
400  transpose matrix, see BuildTranspose().
401 
402  For the serial CPU backends, the internal transpose is not required, and
403  this function is a no-op. This allows for significant memory savings
404  when the internal transpose matrix is not required. */
405  void EnsureMultTranspose() const;
406 
407  void PartMult(const Array<int> &rows, const Vector &x, Vector &y) const;
408  void PartAddMult(const Array<int> &rows, const Vector &x, Vector &y,
409  const double a=1.0) const;
410 
411  /// y = A * x, treating all entries as booleans (zero=false, nonzero=true).
412  /** The actual values stored in the data array, #A, are not used - this means
413  that all entries in the sparsity pattern are considered to be true by
414  this method. */
415  void BooleanMult(const Array<int> &x, Array<int> &y) const;
416 
417  /// y = At * x, treating all entries as booleans (zero=false, nonzero=true).
418  /** The actual values stored in the data array, #A, are not used - this means
419  that all entries in the sparsity pattern are considered to be true by
420  this method. */
421  void BooleanMultTranspose(const Array<int> &x, Array<int> &y) const;
422 
423  /// y = |A| * x, using entry-wise absolute values of matrix A
424  void AbsMult(const Vector &x, Vector &y) const;
425 
426  /// y = |At| * x, using entry-wise absolute values of the transpose of matrix A
427  /** If the matrix is modified, call ResetTranspose() and optionally
428  EnsureMultTranspose() to make sure this method uses the correct updated
429  transpose. */
430  void AbsMultTranspose(const Vector &x, Vector &y) const;
431 
432  /// Compute y^t A x
433  double InnerProduct(const Vector &x, const Vector &y) const;
434 
435  /// For all i compute \f$ x_i = \sum_j A_{ij} \f$
436  void GetRowSums(Vector &x) const;
437  /// For i = irow compute \f$ x_i = \sum_j | A_{i, j} | \f$
438  double GetRowNorml1(int irow) const;
439 
440  /// This virtual method is not supported: it always returns NULL.
441  virtual MatrixInverse *Inverse() const;
442 
443  /// Eliminates a column from the transpose matrix.
444  void EliminateRow(int row, const double sol, Vector &rhs);
445 
446  /// Eliminates a row from the matrix.
447  /*!
448  * - If @a dpolicy = #DIAG_ZERO, all the entries in the row will be set to 0.
449  * - If @a dpolicy = #DIAG_ONE (matrix must be square), the diagonal entry
450  * will be set equal to 1 and all other entries in the row to 0.
451  * - The policy #DIAG_KEEP is not supported.
452  */
453  void EliminateRow(int row, DiagonalPolicy dpolicy = DIAG_ZERO);
454 
455  /// Eliminates the column @a col from the matrix.
456  /** - If @a dpolicy = #DIAG_ZERO, all entries in the column will be set to 0.
457  - If @a dpolicy = #DIAG_ONE (matrix must be square), the diagonal entry
458  will be set equal to 1 and all other entries in the column to 0.
459  - The policy #DIAG_KEEP is not supported. */
460  void EliminateCol(int col, DiagonalPolicy dpolicy = DIAG_ZERO);
461 
462  /// Eliminate all columns i for which @a cols[i] != 0.
463  /** Elimination of a column means that all entries in the column are set to
464  zero. In addition, if the pointers @a x and @a b are not NULL, the
465  eliminated matrix entries are multiplied by the corresponding solution
466  value in @a *x and subtracted from the r.h.s. vector, @a *b. */
467  void EliminateCols(const Array<int> &cols, const Vector *x = NULL,
468  Vector *b = NULL);
469 
470  /** @brief Similar to EliminateCols + save the eliminated entries into
471  @a Ae so that (*this) + Ae is equal to the original matrix. */
472  void EliminateCols(const Array<int> &col_marker, SparseMatrix &Ae);
473 
474  /// Eliminate row @a rc and column @a rc and modify the @a rhs using @a sol.
475  /** Eliminates the column @a rc to the @a rhs, deletes the row @a rc and
476  replaces the element (rc,rc) with 1.0; assumes that element (i,rc)
477  is assembled if and only if the element (rc,i) is assembled.
478  By default, elements (rc,rc) are set to 1.0, although this behavior
479  can be adjusted by changing the @a dpolicy parameter. */
480  void EliminateRowCol(int rc, const double sol, Vector &rhs,
481  DiagonalPolicy dpolicy = DIAG_ONE);
482 
483  /** @brief Similar to
484  EliminateRowCol(int, const double, Vector &, DiagonalPolicy), but
485  multiple values for eliminated unknowns are accepted, and accordingly
486  multiple right-hand-sides are used. */
487  void EliminateRowColMultipleRHS(int rc, const Vector &sol,
488  DenseMatrix &rhs,
489  DiagonalPolicy dpolicy = DIAG_ONE);
490 
491  /// Perform elimination and set the diagonal entry to the given value
492  void EliminateRowColDiag(int rc, double value);
493 
494  /// Eliminate row @a rc and column @a rc.
495  void EliminateRowCol(int rc, DiagonalPolicy dpolicy = DIAG_ONE);
496 
497  /** @brief Similar to EliminateRowCol(int, DiagonalPolicy) + save the
498  eliminated entries into @a Ae so that (*this) + Ae is equal to the
499  original matrix */
500  void EliminateRowCol(int rc, SparseMatrix &Ae,
501  DiagonalPolicy dpolicy = DIAG_ONE);
502 
503  /** @brief Eliminate essential (Dirichlet) boundary conditions.
504 
505  @param[in] ess_dofs indices of the degrees of freedom belonging to the
506  essential boundary conditions.
507  @param[in] diag_policy policy for diagonal entries. */
508  void EliminateBC(const Array<int> &ess_dofs,
509  DiagonalPolicy diag_policy);
510 
511  /// If a row contains only one diag entry of zero, set it to 1.
512  void SetDiagIdentity();
513  /// If a row contains only zeros, set its diagonal to 1.
514  virtual void EliminateZeroRows(const double threshold = 1e-12);
515 
516  /// Gauss-Seidel forward and backward iterations over a vector x.
517  void Gauss_Seidel_forw(const Vector &x, Vector &y) const;
518  void Gauss_Seidel_back(const Vector &x, Vector &y) const;
519 
520  /// Determine appropriate scaling for Jacobi iteration
521  double GetJacobiScaling() const;
522  /** One scaled Jacobi iteration for the system A x = b.
523  x1 = x0 + sc D^{-1} (b - A x0) where D is the diag of A.
524  Absolute values of D are used when use_abs_diag = true. */
525  void Jacobi(const Vector &b, const Vector &x0, Vector &x1,
526  double sc, bool use_abs_diag = false) const;
527 
528  /// x = sc b / A_ii. When use_abs_diag = true, |A_ii| is used.
529  void DiagScale(const Vector &b, Vector &x,
530  double sc = 1.0, bool use_abs_diag = false) const;
531 
532  /** x1 = x0 + sc D^{-1} (b - A x0) where \f$ D_{ii} = \sum_j |A_{ij}| \f$. */
533  void Jacobi2(const Vector &b, const Vector &x0, Vector &x1,
534  double sc = 1.0) const;
535 
536  /** x1 = x0 + sc D^{-1} (b - A x0) where \f$ D_{ii} = \sum_j A_{ij} \f$. */
537  void Jacobi3(const Vector &b, const Vector &x0, Vector &x1,
538  double sc = 1.0) const;
539 
540  /** @brief Finalize the matrix initialization, switching the storage format
541  from LIL to CSR. */
542  /** This method should be called once, after the matrix has been initialized.
543  Internally, this method converts the matrix from row-wise linked list
544  (LIL) format into CSR (compressed sparse row) format. */
545  virtual void Finalize(int skip_zeros = 1) { Finalize(skip_zeros, false); }
546 
547  /// A slightly more general version of the Finalize(int) method.
548  void Finalize(int skip_zeros, bool fix_empty_rows);
549 
550  /// Returns whether or not CSR format has been finalized.
551  bool Finalized() const { return !A.Empty(); }
552  /// Returns whether or not the columns are sorted.
553  bool ColumnsAreSorted() const { return isSorted; }
554 
555  /** @brief Remove entries smaller in absolute value than a given tolerance
556  @a tol. If @a fix_empty_rows is true, a zero value is inserted in the
557  diagonal entry (for square matrices only) */
558  void Threshold(double tol, bool fix_empty_rows = false);
559 
560  /** Split the matrix into M x N blocks of sparse matrices in CSR format.
561  The 'blocks' array is M x N (i.e. M and N are determined by its
562  dimensions) and its entries are overwritten by the new blocks. */
563  void GetBlocks(Array2D<SparseMatrix *> &blocks) const;
564 
565  void GetSubMatrix(const Array<int> &rows, const Array<int> &cols,
566  DenseMatrix &subm) const;
567 
568  /** @brief Initialize the SparseMatrix for fast access to the entries of the
569  given @a row which becomes the "current row". */
570  /** Fast access to the entries of the "current row" can be performed using
571  the methods: SearchRow(const int), _Add_(const int, const double),
572  _Set_(const int, const double), and _Get_(const int). */
573  inline void SetColPtr(const int row) const;
574  /** @brief Reset the "current row" set by calling SetColPtr(). This method
575  must be called between any two calls to SetColPtr(). */
576  inline void ClearColPtr() const;
577  /// Perform a fast search for an entry in the "current row". See SetColPtr().
578  /** If the matrix is not finalized and the entry is not found in the
579  SparseMatrix, it will be added to the sparsity pattern initialized with
580  zero. If the matrix is finalized and the entry is not found, an error
581  will be generated. */
582  inline double &SearchRow(const int col);
583  /// Add a value to an entry in the "current row". See SetColPtr().
584  inline void _Add_(const int col, const double a)
585  { SearchRow(col) += a; }
586  /// Set an entry in the "current row". See SetColPtr().
587  inline void _Set_(const int col, const double a)
588  { SearchRow(col) = a; }
589  /// Read the value of an entry in the "current row". See SetColPtr().
590  inline double _Get_(const int col) const;
591 
592  inline double &SearchRow(const int row, const int col);
593  inline void _Add_(const int row, const int col, const double a)
594  { SearchRow(row, col) += a; }
595  inline void _Set_(const int row, const int col, const double a)
596  { SearchRow(row, col) = a; }
597 
598  void Set(const int i, const int j, const double val);
599  void Add(const int i, const int j, const double val);
600 
601  void SetSubMatrix(const Array<int> &rows, const Array<int> &cols,
602  const DenseMatrix &subm, int skip_zeros = 1);
603 
604  void SetSubMatrixTranspose(const Array<int> &rows, const Array<int> &cols,
605  const DenseMatrix &subm, int skip_zeros = 1);
606 
607  /** Insert the DenseMatrix into this SparseMatrix at the specified rows and
608  columns. If \c skip_zeros==0 , all entries from the DenseMatrix are
609  added including zeros. If \c skip_zeros==2 , no zeros are added to the
610  SparseMatrix regardless of their position in the matrix. Otherwise, the
611  default \c skip_zeros behavior is to omit the zero from the SparseMatrix
612  unless it would break the symmetric structure of the SparseMatrix. */
613  void AddSubMatrix(const Array<int> &rows, const Array<int> &cols,
614  const DenseMatrix &subm, int skip_zeros = 1);
615 
616  bool RowIsEmpty(const int row) const;
617 
618  /// Extract all column indices and values from a given row.
619  /** If the matrix is finalized (i.e. in CSR format), @a cols and @a srow will
620  simply be references to the specific portion of the #J and #A arrays.
621  As required by the AbstractSparseMatrix interface this method returns:
622  - 0, if @a cols and @a srow are copies of the values in the matrix, i.e.
623  when the matrix is open.
624  - 1, if @a cols and @a srow are views of the values in the matrix, i.e.
625  when the matrix is finalized.
626  @warning This method breaks the const-ness when the matrix is finalized
627  because it gives write access to the #J and #A arrays. */
628  virtual int GetRow(const int row, Array<int> &cols, Vector &srow) const;
629 
630  void SetRow(const int row, const Array<int> &cols, const Vector &srow);
631  void AddRow(const int row, const Array<int> &cols, const Vector &srow);
632 
633  void ScaleRow(const int row, const double scale);
634  /// this = diag(sl) * this;
635  void ScaleRows(const Vector & sl);
636  /// this = this * diag(sr);
637  void ScaleColumns(const Vector & sr);
638 
639  /** @brief Add the sparse matrix 'B' to '*this'. This operation will cause an
640  error if '*this' is finalized and 'B' has larger sparsity pattern. */
642 
643  /** @brief Add the sparse matrix 'B' scaled by the scalar 'a' into '*this'.
644  Only entries in the sparsity pattern of '*this' are added. */
645  void Add(const double a, const SparseMatrix &B);
646 
647  SparseMatrix &operator=(double a);
648 
649  SparseMatrix &operator*=(double a);
650 
651  /// Prints matrix to stream out.
652  void Print(std::ostream &out = mfem::out, int width_ = 4) const;
653 
654  /// Prints matrix in matlab format.
655  virtual void PrintMatlab(std::ostream &out = mfem::out) const;
656 
657  /// Prints matrix in Matrix Market sparse format.
658  void PrintMM(std::ostream &out = mfem::out) const;
659 
660  /// Prints matrix to stream out in hypre_CSRMatrix format.
661  void PrintCSR(std::ostream &out) const;
662 
663  /// Prints a sparse matrix to stream out in CSR format.
664  void PrintCSR2(std::ostream &out) const;
665 
666  /// Print various sparse matrix statistics.
667  void PrintInfo(std::ostream &out) const;
668 
669  /// Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix
670  double IsSymmetric() const;
671 
672  /// (*this) = 1/2 ((*this) + (*this)^t)
673  void Symmetrize();
674 
675  /// Returns the number of the nonzero elements in the matrix
676  virtual int NumNonZeroElems() const;
677 
678  double MaxNorm() const;
679 
680  /// Count the number of entries with |a_ij| <= tol.
681  int CountSmallElems(double tol) const;
682 
683  /// Count the number of entries that are NOT finite, i.e. Inf or Nan.
684  int CheckFinite() const;
685 
686  /// Set the graph ownership flag (I and J arrays).
687  void SetGraphOwner(bool ownij)
688  { I.SetHostPtrOwner(ownij); J.SetHostPtrOwner(ownij); }
689 
690  /// Set the data ownership flag (A array).
691  void SetDataOwner(bool owna) { A.SetHostPtrOwner(owna); }
692 
693  /// Get the graph ownership flag (I and J arrays).
694  bool OwnsGraph() const { return I.OwnsHostPtr() && J.OwnsHostPtr(); }
695 
696  /// Get the data ownership flag (A array).
697  bool OwnsData() const { return A.OwnsHostPtr(); }
698 
699  /// Lose the ownership of the graph (I, J) and data (A) arrays.
700  void LoseData() { SetGraphOwner(false); SetDataOwner(false); }
701 
702  void Swap(SparseMatrix &other);
703 
704  /// Destroys sparse matrix.
705  virtual ~SparseMatrix();
706 
707  Type GetType() const { return MFEM_SPARSEMAT; }
708 };
709 
710 inline std::ostream& operator<<(std::ostream& os, SparseMatrix const& mat)
711 {
712  mat.Print(os);
713  return os;
714 }
715 
716 /// Applies f() to each element of the matrix (after it is finalized).
717 void SparseMatrixFunction(SparseMatrix &S, double (*f)(double));
718 
719 
720 /// Transpose of a sparse matrix. A must be finalized.
721 SparseMatrix *Transpose(const SparseMatrix &A);
722 /// Transpose of a sparse matrix. A does not need to be a CSR matrix.
723 SparseMatrix *TransposeAbstractSparseMatrix (const AbstractSparseMatrix &A,
724  int useActualWidth);
725 
726 /// Matrix product A.B.
727 /** If @a OAB is not NULL, we assume it has the structure of A.B and store the
728  result in @a OAB. If @a OAB is NULL, we create a new SparseMatrix to store
729  the result and return a pointer to it.
730 
731  All matrices must be finalized. */
732 SparseMatrix *Mult(const SparseMatrix &A, const SparseMatrix &B,
733  SparseMatrix *OAB = NULL);
734 
735 /// C = A^T B
736 SparseMatrix *TransposeMult(const SparseMatrix &A, const SparseMatrix &B);
737 
738 /// Matrix product of sparse matrices. A and B do not need to be CSR matrices
739 SparseMatrix *MultAbstractSparseMatrix (const AbstractSparseMatrix &A,
740  const AbstractSparseMatrix &B);
741 
742 /// Matrix product A.B
743 DenseMatrix *Mult(const SparseMatrix &A, DenseMatrix &B);
744 
745 /// RAP matrix product (with R=P^T)
746 DenseMatrix *RAP(const SparseMatrix &A, DenseMatrix &P);
747 
748 /// RAP matrix product (with R=P^T)
749 DenseMatrix *RAP(DenseMatrix &A, const SparseMatrix &P);
750 
751 /** RAP matrix product (with P=R^T). ORAP is like OAB above.
752  All matrices must be finalized. */
753 SparseMatrix *RAP(const SparseMatrix &A, const SparseMatrix &R,
754  SparseMatrix *ORAP = NULL);
755 
756 /// General RAP with given R^T, A and P
757 SparseMatrix *RAP(const SparseMatrix &Rt, const SparseMatrix &A,
758  const SparseMatrix &P);
759 
760 /// Matrix multiplication A^t D A. All matrices must be finalized.
761 SparseMatrix *Mult_AtDA(const SparseMatrix &A, const Vector &D,
762  SparseMatrix *OAtDA = NULL);
763 
764 
765 /// Matrix addition result = A + B.
766 SparseMatrix * Add(const SparseMatrix & A, const SparseMatrix & B);
767 /// Matrix addition result = a*A + b*B
768 SparseMatrix * Add(double a, const SparseMatrix & A, double b,
769  const SparseMatrix & B);
770 /// Matrix addition result = sum_i A_i
771 SparseMatrix * Add(Array<SparseMatrix *> & Ai);
772 
773 /// B += alpha * A
774 void Add(const SparseMatrix &A, double alpha, DenseMatrix &B);
775 
776 /// Produces a block matrix with blocks A_{ij}*B
777 DenseMatrix *OuterProduct(const DenseMatrix &A, const DenseMatrix &B);
778 
779 /// Produces a block matrix with blocks A_{ij}*B
780 SparseMatrix *OuterProduct(const DenseMatrix &A, const SparseMatrix &B);
781 
782 /// Produces a block matrix with blocks A_{ij}*B
783 SparseMatrix *OuterProduct(const SparseMatrix &A, const DenseMatrix &B);
784 
785 /// Produces a block matrix with blocks A_{ij}*B
786 SparseMatrix *OuterProduct(const SparseMatrix &A, const SparseMatrix &B);
787 
788 
789 // Inline methods
790 
791 inline void SparseMatrix::SetColPtr(const int row) const
792 {
793  if (Rows)
794  {
795  if (ColPtrNode == NULL)
796  {
797  ColPtrNode = new RowNode *[width];
798  for (int i = 0; i < width; i++)
799  {
800  ColPtrNode[i] = NULL;
801  }
802  }
803  for (RowNode *node_p = Rows[row]; node_p != NULL; node_p = node_p->Prev)
804  {
805  ColPtrNode[node_p->Column] = node_p;
806  }
807  }
808  else
809  {
810  if (ColPtrJ == NULL)
811  {
812  ColPtrJ = new int[width];
813  for (int i = 0; i < width; i++)
814  {
815  ColPtrJ[i] = -1;
816  }
817  }
818  for (int j = I[row], end = I[row+1]; j < end; j++)
819  {
820  ColPtrJ[J[j]] = j;
821  }
822  }
823  current_row = row;
824 }
825 
826 inline void SparseMatrix::ClearColPtr() const
827 {
828  if (Rows)
829  {
830  for (RowNode *node_p = Rows[current_row]; node_p != NULL;
831  node_p = node_p->Prev)
832  {
833  ColPtrNode[node_p->Column] = NULL;
834  }
835  }
836  else
837  {
838  for (int j = I[current_row], end = I[current_row+1]; j < end; j++)
839  {
840  ColPtrJ[J[j]] = -1;
841  }
842  }
843 }
844 
845 inline double &SparseMatrix::SearchRow(const int col)
846 {
847  if (Rows)
848  {
849  RowNode *node_p = ColPtrNode[col];
850  if (node_p == NULL)
851  {
852 #ifdef MFEM_USE_MEMALLOC
853  node_p = NodesMem->Alloc();
854 #else
855  node_p = new RowNode;
856 #endif
857  node_p->Prev = Rows[current_row];
858  node_p->Column = col;
859  node_p->Value = 0.0;
860  Rows[current_row] = ColPtrNode[col] = node_p;
861  }
862  return node_p->Value;
863  }
864  else
865  {
866  const int j = ColPtrJ[col];
867  MFEM_VERIFY(j != -1, "Entry for column " << col << " is not allocated.");
868  return A[j];
869  }
870 }
871 
872 inline double SparseMatrix::_Get_(const int col) const
873 {
874  if (Rows)
875  {
876  RowNode *node_p = ColPtrNode[col];
877  return (node_p == NULL) ? 0.0 : node_p->Value;
878  }
879  else
880  {
881  const int j = ColPtrJ[col];
882  return (j == -1) ? 0.0 : A[j];
883  }
884 }
885 
886 inline double &SparseMatrix::SearchRow(const int row, const int col)
887 {
888  if (Rows)
889  {
890  RowNode *node_p;
891 
892  for (node_p = Rows[row]; 1; node_p = node_p->Prev)
893  {
894  if (node_p == NULL)
895  {
896 #ifdef MFEM_USE_MEMALLOC
897  node_p = NodesMem->Alloc();
898 #else
899  node_p = new RowNode;
900 #endif
901  node_p->Prev = Rows[row];
902  node_p->Column = col;
903  node_p->Value = 0.0;
904  Rows[row] = node_p;
905  break;
906  }
907  else if (node_p->Column == col)
908  {
909  break;
910  }
911  }
912  return node_p->Value;
913  }
914  else
915  {
916  int *Ip = I+row, *Jp = J;
917  for (int k = Ip[0], end = Ip[1]; k < end; k++)
918  {
919  if (Jp[k] == col)
920  {
921  return A[k];
922  }
923  }
924  MFEM_ABORT("Could not find entry for row = " << row << ", col = " << col);
925  }
926  return A[0];
927 }
928 
929 /// Specialization of the template function Swap<> for class SparseMatrix
930 template<> inline void Swap<SparseMatrix>(SparseMatrix &a, SparseMatrix &b)
931 {
932  a.Swap(b);
933 }
934 
935 } // namespace mfem
936 
937 #endif
RowNode ** ColPtrNode
Definition: sparsemat.hpp:77
Memory< int > I
Array with size (height+1) containing the row offsets.
Definition: sparsemat.hpp:62
Elem * Alloc()
Definition: mem_alloc.hpp:166
void _Add_(const int row, const int col, const double a)
Definition: sparsemat.hpp:593
double GetJacobiScaling() const
Determine appropriate scaling for Jacobi iteration.
Definition: sparsemat.cpp:2551
int CheckFinite() const
Count the number of entries that are NOT finite, i.e. Inf or Nan.
Definition: sparsemat.cpp:1660
cusparseStatus_t status
Definition: sparsemat.hpp:106
void DiagScale(const Vector &b, Vector &x, double sc=1.0, bool use_abs_diag=false) const
x = sc b / A_ii. When use_abs_diag = true, |A_ii| is used.
Definition: sparsemat.cpp:2616
SparseMatrix * TransposeAbstractSparseMatrix(const AbstractSparseMatrix &A, int useActualWidth)
Transpose of a sparse matrix. A does not need to be a CSR matrix.
Definition: sparsemat.cpp:3558
int RowSize(const int i) const
Returns the number of elements in row i.
Definition: sparsemat.cpp:344
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:1584
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:584
void BuildTranspose() const
Build and store internally the transpose of this matrix which will be used in the methods AddMultTran...
Definition: sparsemat.cpp:956
void Clear()
Clear the contents of the SparseMatrix.
Definition: sparsemat.hpp:209
void EliminateCol(int col, DiagonalPolicy dpolicy=DIAG_ZERO)
Eliminates the column col from the matrix.
Definition: sparsemat.cpp:1737
SparseMatrix * At
Transpose of A. Owned. Used to perform MultTranspose() on devices.
Definition: sparsemat.hpp:80
Type GetType() const
Definition: sparsemat.hpp:707
static MemoryClass GetHostMemoryClass()
Get the current Host MemoryClass. This is the MemoryClass used by most MFEM host Memory objects...
Definition: device.hpp:268
bool ColumnsAreSorted() const
Returns whether or not the columns are sorted.
Definition: sparsemat.hpp:553
void _Set_(const int row, const int col, const double a)
Definition: sparsemat.hpp:595
static size_t bufferSize
Definition: sparsemat.hpp:101
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
Definition: sparsemat.hpp:545
cusparseMatDescr_t descr
Definition: sparsemat.hpp:108
void MakeRef(const SparseMatrix &master)
Clear the contents of the SparseMatrix and make it a reference to master.
Definition: sparsemat.cpp:313
double & SearchRow(const int col)
Perform a fast search for an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:845
DenseMatrix * ToDenseMatrix() const
Produces a DenseMatrix from a SparseMatrix.
Definition: sparsemat.cpp:699
bool Empty() const
Check if the SparseMatrix is empty.
Definition: sparsemat.hpp:219
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Definition: sparsemat.cpp:391
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:1022
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
const double * ReadData(bool on_dev=true) const
Definition: sparsemat.hpp:271
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:227
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:3137
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:336
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:736
double MaxNorm() const
Definition: sparsemat.cpp:1607
int * GetI()
Return the array I.
Definition: sparsemat.hpp:222
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:911
Abstract data type for matrix inverse.
Definition: matrix.hpp:62
static void * dBuffer
Definition: sparsemat.hpp:102
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:1949
void ClearGPUSparse()
Clear the cuSPARSE/hipSPARSE descriptors. This must be called after releasing the device memory of A...
Definition: sparsemat.cpp:81
void PrintInfo(std::ostream &out) const
Print various sparse matrix statistics.
Definition: sparsemat.cpp:3381
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
Definition: sparsemat.cpp:597
void Jacobi(const Vector &b, const Vector &x0, Vector &x1, double sc, bool use_abs_diag=false) const
Definition: sparsemat.cpp:2584
void LoseData()
Lose the ownership of the graph (I, J) and data (A) arrays.
Definition: sparsemat.hpp:700
void AddRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:3018
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
Definition: sparsemat.cpp:405
void SortColumnIndices()
Sort the column indices corresponding to each row.
Definition: sparsemat.cpp:457
Memory< double > & GetMemoryData()
Definition: sparsemat.hpp:269
void PrintMM(std::ostream &out=mfem::out) const
Prints matrix in Matrix Market sparse format.
Definition: sparsemat.cpp:3310
int Capacity() const
Return the size of the allocated memory.
Memory< int > & GetMemoryI()
Definition: sparsemat.hpp:237
cusparseDnVecDescr_t vecY_descr
Definition: sparsemat.hpp:113
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const
Definition: sparsemat.cpp:2882
double * HostReadWriteData()
Definition: sparsemat.hpp:281
void _Set_(const int col, const double a)
Set an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:587
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:1052
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:3926
void PartAddMult(const Array< int > &rows, const Vector &x, Vector &y, const double a=1.0) const
Definition: sparsemat.cpp:1004
bool RowIsEmpty(const int row) const
Definition: sparsemat.cpp:2910
void GetBlocks(Array2D< SparseMatrix * > &blocks) const
Definition: sparsemat.cpp:1433
const int * ReadJ(bool on_dev=true) const
Definition: sparsemat.hpp:255
void ScaleRow(const int row, const double scale)
Definition: sparsemat.cpp:3050
RowNodeAlloc * NodesMem
Definition: sparsemat.hpp:84
double * GetData()
Return the element data, i.e. the array A.
Definition: sparsemat.hpp:232
Memory< int > & GetMemoryJ()
Definition: sparsemat.hpp:253
void UseGPUSparse(bool useGPUSparse_=true)
Runtime option to use cuSPARSE or hipSPARSE. Only valid when using a CUDA or HIP backend.
Definition: sparsemat.hpp:191
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
Definition: sparsemat.cpp:1565
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:573
void SetGraphOwner(bool ownij)
Set the graph ownership flag (I and J arrays).
Definition: sparsemat.hpp:687
ID for class SparseMatrix.
Definition: operator.hpp:259
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2282
void Set(const int i, const int j, const double val)
Definition: sparsemat.cpp:2764
hipsparseDnVecDescr_t vecY_descr
Definition: sparsemat.hpp:125
bool OwnsData() const
Get the data ownership flag (A array).
Definition: sparsemat.hpp:697
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:826
MFEM_DEPRECATED void ClearCuSparse()
Deprecated equivalent of ClearGPUSparse().
Definition: sparsemat.hpp:216
void ScaleColumns(const Vector &sr)
this = this * diag(sr);
Definition: sparsemat.cpp:3109
virtual ~SparseMatrix()
Destroys sparse matrix.
Definition: sparsemat.cpp:4198
void PrintCSR2(std::ostream &out) const
Prints a sparse matrix to stream out in CSR format.
Definition: sparsemat.cpp:3356
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
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:2931
bool isSorted
Are the columns sorted already.
Definition: sparsemat.hpp:88
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:68
SparseMatrix()
Create an empty SparseMatrix.
Definition: sparsemat.hpp:131
Data type sparse matrix.
Definition: sparsemat.hpp:50
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
void EnsureMultTranspose() const
Ensures that the matrix is capable of performing MultTranspose(), AddMultTranspose(), and AbsMultTranspose().
Definition: sparsemat.cpp:970
static MemoryClass GetDeviceMemoryClass()
Get the current Device MemoryClass. This is the MemoryClass used by most MFEM device kernels to acces...
Definition: device.hpp:281
hipsparseStatus_t status
Definition: sparsemat.hpp:119
hipsparseSpMatDescr_t matA_descr
Definition: sparsemat.hpp:123
bool OwnsGraph() const
Get the graph ownership flag (I and J arrays).
Definition: sparsemat.hpp:694
void SetWidth(int width_=-1)
Change the width of a SparseMatrix.
Definition: sparsemat.cpp:419
double b
Definition: lissajous.cpp:42
const double * HostReadData() const
Definition: sparsemat.hpp:277
int * ReadWriteI(bool on_dev=true)
Definition: sparsemat.hpp:243
double IsSymmetric() const
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix.
Definition: sparsemat.cpp:1523
const int * HostReadJ() const
Definition: sparsemat.hpp:261
static int SparseMatrixCount
Definition: sparsemat.hpp:100
void SetDiagIdentity()
If a row contains only one diag entry of zero, set it to 1.
Definition: sparsemat.cpp:2354
MFEM_DEPRECATED void UseCuSparse(bool useCuSparse_=true)
Deprecated equivalent of UseGPUSparse().
Definition: sparsemat.hpp:194
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:3078
static cusparseHandle_t handle
Definition: sparsemat.hpp:107
void Swap< SparseMatrix >(SparseMatrix &a, SparseMatrix &b)
Specialization of the template function Swap&lt;&gt; for class SparseMatrix.
Definition: sparsemat.hpp:930
DenseMatrix * OuterProduct(const DenseMatrix &A, const DenseMatrix &B)
Produces a block matrix with blocks A_{ij}*B.
Definition: sparsemat.cpp:4070
int * WriteI(bool on_dev=true)
Definition: sparsemat.hpp:241
Set the diagonal value to one.
Definition: operator.hpp:50
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
Definition: densemat.cpp:3213
virtual void PrintMatlab(std::ostream &out=mfem::out) const
Prints matrix in matlab format.
Definition: sparsemat.cpp:3288
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:319
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2718
void Add(const int i, const int j, const double val)
Definition: sparsemat.cpp:2783
Dynamic 2D array using row-major layout.
Definition: array.hpp:351
void SetSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2802
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
Definition: sparsemat.hpp:342
bool Finalized() const
Returns whether or not CSR format has been finalized.
Definition: sparsemat.hpp:551
const double * GetData() const
Return the element data, i.e. the array A, const version.
Definition: sparsemat.hpp:234
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:256
void SparseMatrixFunction(SparseMatrix &S, double(*f)(double))
Applies f() to each element of the matrix (after it is finalized).
Definition: sparsemat.cpp:3491
virtual MatrixInverse * Inverse() const
This virtual method is not supported: it always returns NULL.
Definition: sparsemat.cpp:1684
const Memory< double > & GetMemoryData() const
Definition: sparsemat.hpp:270
void EliminateRowColDiag(int rc, double value)
Perform elimination and set the diagonal entry to the given value.
Definition: sparsemat.cpp:2154
const Memory< int > & GetMemoryI() const
Definition: sparsemat.hpp:238
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:729
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:791
const int * GetJ() const
Return the array J, const version.
Definition: sparsemat.hpp:229
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
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:2841
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:1850
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Definition: table.cpp:413
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:1257
double * HostWriteData()
Definition: sparsemat.hpp:279
int Size() const
For backward compatibility, define Size() to be synonym of Height().
Definition: sparsemat.hpp:206
void Gauss_Seidel_back(const Vector &x, Vector &y) const
Definition: sparsemat.cpp:2468
const int * ReadI(bool on_dev=true) const
Definition: sparsemat.hpp:239
void Gauss_Seidel_forw(const Vector &x, Vector &y) const
Gauss-Seidel forward and backward iterations over a vector x.
Definition: sparsemat.cpp:2384
int * HostReadWriteJ()
Definition: sparsemat.hpp:265
RowNode ** Rows
Array of linked lists, one for every row. This array represents the linked list (LIL) storage format...
Definition: sparsemat.hpp:73
hipsparseDnVecDescr_t vecX_descr
Definition: sparsemat.hpp:124
const Memory< int > & GetMemoryJ() const
Definition: sparsemat.hpp:254
int * WriteJ(bool on_dev=true)
Definition: sparsemat.hpp:257
void EliminateRow(int row, const double sol, Vector &rhs)
Eliminates a column from the transpose matrix.
Definition: sparsemat.cpp:1689
double a
Definition: lissajous.cpp:41
int ActualWidth() const
Returns the actual Width of the matrix.
Definition: sparsemat.cpp:3465
void GetDiag(Vector &d) const
Returns the Diagonal of A.
Definition: sparsemat.cpp:665
MemAlloc< RowNode, 1024 > RowNodeAlloc
Definition: sparsemat.hpp:83
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:353
double GetRowNorml1(int irow) const
For i = irow compute .
Definition: sparsemat.cpp:1233
void ResetTranspose() const
Definition: sparsemat.cpp:964
SparseMatrix * TransposeMult(const SparseMatrix &A, const SparseMatrix &B)
C = A^T B.
Definition: sparsemat.cpp:3763
virtual void EliminateZeroRows(const double threshold=1e-12)
If a row contains only zeros, set its diagonal to 1.
Definition: sparsemat.cpp:2365
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:1124
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
Definition: operator.hpp:47
double & operator()(int i, int j)
Returns reference to A[i][j].
Definition: sparsemat.cpp:607
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:904
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:1777
void OverrideSize(int height_, int width_)
Sets the height and width of the matrix.
Definition: sparsemat.cpp:297
SparseMatrix & operator=(const SparseMatrix &rhs)
Assignment operator: deep copy.
Definition: sparsemat.cpp:303
SparseMatrix & operator*=(double a)
Definition: sparsemat.cpp:3216
void SetDataOwner(bool owna)
Set the data ownership flag (A array).
Definition: sparsemat.hpp:691
std::ostream & operator<<(std::ostream &os, SparseMatrix const &mat)
Definition: sparsemat.hpp:710
void SetRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:2971
cusparseSpMatDescr_t matA_descr
Definition: sparsemat.hpp:111
const double alpha
Definition: ex15.cpp:369
double * ReadWriteData(bool on_dev=true)
Definition: sparsemat.hpp:275
int * ReadWriteJ(bool on_dev=true)
Definition: sparsemat.hpp:259
static hipsparseHandle_t handle
Definition: sparsemat.hpp:120
Vector data type.
Definition: vector.hpp:60
void Print(std::ostream &out=mfem::out, int width_=4) const
Prints matrix to stream out.
Definition: sparsemat.cpp:3240
const int * HostReadI() const
Definition: sparsemat.hpp:245
const int * GetI() const
Return the array I, const version.
Definition: sparsemat.hpp:224
void Jacobi2(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
Definition: sparsemat.cpp:2704
cusparseDnVecDescr_t vecX_descr
Definition: sparsemat.hpp:112
int MaxRowSize() const
Returns the maximum number of elements among all rows.
Definition: sparsemat.cpp:367
int * HostReadWriteI()
Definition: sparsemat.hpp:249
void Jacobi3(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
Definition: sparsemat.cpp:2711
cusparseMatDescr_t matA_descr
Definition: sparsemat.hpp:115
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:3771
void AbsMult(const Vector &x, Vector &y) const
y = |A| * x, using entry-wise absolute values of matrix A
Definition: sparsemat.cpp:1075
double _Get_(const int col) const
Read the value of an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:872
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:65
void Swap(SparseMatrix &other)
Definition: sparsemat.cpp:4178
void PrintCSR(std::ostream &out) const
Prints matrix to stream out in hypre_CSRMatrix format.
Definition: sparsemat.cpp:3332
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:1169
MemoryClass
Memory classes identify sets of memory types.
Definition: mem_manager.hpp:73
Set the diagonal value to zero.
Definition: operator.hpp:49
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
void EliminateBC(const Array< int > &ess_dofs, DiagonalPolicy diag_policy)
Eliminate essential (Dirichlet) boundary conditions.
Definition: sparsemat.cpp:2311
int CountSmallElems(double tol) const
Count the number of entries with |a_ij| &lt;= tol.
Definition: sparsemat.cpp:1632
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const
Definition: sparsemat.cpp:978
void GetRowSums(Vector &x) const
For all i compute .
Definition: sparsemat.cpp:1210
double * WriteData(bool on_dev=true)
Definition: sparsemat.hpp:273