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