MFEM  v4.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, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_SPARSEMAT_HPP
13 #define MFEM_SPARSEMAT_HPP
14 
15 // Data types for sparse matrix
16 
17 #include "../general/mem_alloc.hpp"
18 #include "../general/mem_manager.hpp"
19 #include "../general/device.hpp"
20 #include "../general/table.hpp"
21 #include "../general/globals.hpp"
22 #include "densemat.hpp"
23 
24 namespace mfem
25 {
26 
27 class
28 #if defined(__alignas_is_defined)
29  alignas(double)
30 #endif
31  RowNode
32 {
33 public:
34  double Value;
35  RowNode *Prev;
36  int Column;
37 };
38 
39 /// Data type sparse matrix
41 {
42 protected:
43  /// @name Arrays used by the CSR storage format.
44  /** */
45  ///@{
46  /// @brief %Array with size (#height+1) containing the row offsets.
47  /** The data for row r, 0 <= r < height, is at offsets j, I[r] <= j < I[r+1].
48  The offsets, j, are indices in the #J and #A arrays. The first entry in
49  this array is always zero, I[0] = 0, and the last entry, I[height], gives
50  the total number of entries stored (at a minimum, all nonzeros must be
51  represented) in the sparse matrix. */
53  /** @brief %Array with size #I[#height], containing the column indices for
54  all matrix entries, as indexed by the #I array. */
56  /** @brief %Array with size #I[#height], containing the actual entries of the
57  sparse matrix, as indexed by the #I array. */
59  ///@}
60 
61  /** @brief %Array of linked lists, one for every row. This array represents
62  the linked list (LIL) storage format. */
63  RowNode **Rows;
64 
65  mutable int current_row;
66  mutable int* ColPtrJ;
67  mutable RowNode ** ColPtrNode;
68 
69  /// Transpose of A. Owned. Used to perform MultTranspose() on devices.
70  mutable SparseMatrix *At;
71 
72 #ifdef MFEM_USE_MEMALLOC
75 #endif
76 
77  /// Are the columns sorted already.
78  bool isSorted;
79 
80  void Destroy(); // Delete all owned data
81  void SetEmpty(); // Init all entries with empty values
82 
83 public:
84  /// Create an empty SparseMatrix.
86 
87  /** @brief Create a sparse matrix with flexible sparsity structure using a
88  row-wise linked list (LIL) format. */
89  /** New entries are added as needed by methods like AddSubMatrix(),
90  SetSubMatrix(), etc. Calling Finalize() will convert the SparseMatrix to
91  the more compact compressed sparse row (CSR) format. */
92  explicit SparseMatrix(int nrows, int ncols = -1);
93 
94  /** @brief Create a sparse matrix in CSR format. Ownership of @a i, @a j, and
95  @a data is transferred to the SparseMatrix. */
96  SparseMatrix(int *i, int *j, double *data, int m, int n);
97 
98  /** @brief Create a sparse matrix in CSR format. Ownership of @a i, @a j, and
99  @a data is optionally transferred to the SparseMatrix. */
100  /** If the parameter @a data is NULL, then the internal #A array is allocated
101  by this constructor (initializing it with zeros and taking ownership,
102  regardless of the parameter @a owna). */
103  SparseMatrix(int *i, int *j, double *data, int m, int n, bool ownij,
104  bool owna, bool issorted);
105 
106  /** @brief Create a sparse matrix in CSR format where each row has space
107  allocated for exactly @a rowsize entries. */
108  /** SetRow() can then be called or the #I, #J, #A arrays can be used
109  directly. */
110  SparseMatrix(int nrows, int ncols, int rowsize);
111 
112  /// Copy constructor (deep copy).
113  /** If @a mat is finalized and @a copy_graph is false, the #I and #J arrays
114  will use a shallow copy (copy the pointers only) without transferring
115  ownership. */
116  SparseMatrix(const SparseMatrix &mat, bool copy_graph = true);
117 
118  /// Create a SparseMatrix with diagonal @a v, i.e. A = Diag(v)
119  SparseMatrix(const Vector & v);
120 
121 
122  /// Assignment operator: deep copy
123  SparseMatrix& operator=(const SparseMatrix &rhs);
124 
125  /** @brief Clear the contents of the SparseMatrix and make it a reference to
126  @a master */
127  /** After this call, the matrix will point to the same data as @a master but
128  it will not own its data. The @a master must be finalized. */
129  void MakeRef(const SparseMatrix &master);
130 
131  /// For backward compatibility, define Size() to be synonym of Height().
132  int Size() const { return Height(); }
133 
134  /// Clear the contents of the SparseMatrix.
135  void Clear() { Destroy(); SetEmpty(); }
136 
137  /// Check if the SparseMatrix is empty.
138  bool Empty() const { return (A == NULL) && (Rows == NULL); }
139 
140  /// Return the array #I.
141  inline int *GetI() { return I; }
142  /// Return the array #I, const version.
143  inline const int *GetI() const { return I; }
144 
145  /// Return the array #J.
146  inline int *GetJ() { return J; }
147  /// Return the array #J, const version.
148  inline const int *GetJ() const { return J; }
149 
150  /// Return the element data, i.e. the array #A.
151  inline double *GetData() { return A; }
152  /// Return the element data, i.e. the array #A, const version.
153  inline const double *GetData() const { return A; }
154 
155  /// Returns the number of elements in row @a i.
156  int RowSize(const int i) const;
157 
158  /// Returns the maximum number of elements among all rows.
159  int MaxRowSize() const;
160 
161  /// Return a pointer to the column indices in a row.
162  int *GetRowColumns(const int row);
163  /// Return a pointer to the column indices in a row, const version.
164  const int *GetRowColumns(const int row) const;
165 
166  /// Return a pointer to the entries in a row.
167  double *GetRowEntries(const int row);
168  /// Return a pointer to the entries in a row, const version.
169  const double *GetRowEntries(const int row) const;
170 
171  /// Change the width of a SparseMatrix.
172  /*!
173  * If width_ = -1 (DEFAULT), this routine will set the new width
174  * to the actual Width of the matrix awidth = max(J) + 1.
175  * Values 0 <= width_ < awidth are not allowed (error check in Debug Mode only)
176  *
177  * This method can be called for matrices finalized or not.
178  */
179  void SetWidth(int width_ = -1);
180 
181  /// Returns the actual Width of the matrix.
182  /*! This method can be called for matrices finalized or not. */
183  int ActualWidth() const;
184 
185  /// Sort the column indices corresponding to each row.
186  void SortColumnIndices();
187 
188  /** @brief Move the diagonal entry to the first position in each row,
189  preserving the order of the rest of the columns. */
190  void MoveDiagonalFirst();
191 
192  /// Returns reference to a_{ij}.
193  virtual double &Elem(int i, int j);
194 
195  /// Returns constant reference to a_{ij}.
196  virtual const double &Elem(int i, int j) const;
197 
198  /// Returns reference to A[i][j].
199  double &operator()(int i, int j);
200 
201  /// Returns reference to A[i][j].
202  const double &operator()(int i, int j) const;
203 
204  /// Returns the Diagonal of A
205  void GetDiag(Vector & d) const;
206 
207  /// Produces a DenseMatrix from a SparseMatrix
208  DenseMatrix *ToDenseMatrix() const;
209 
210  /// Produces a DenseMatrix from a SparseMatrix
211  void ToDenseMatrix(DenseMatrix & B) const;
212 
213  virtual MemoryClass GetMemoryClass() const
215 
216  /// Matrix vector multiplication.
217  virtual void Mult(const Vector &x, Vector &y) const;
218 
219  /// y += A * x (default) or y += a * A * x
220  void AddMult(const Vector &x, Vector &y, const double a = 1.0) const;
221 
222  /// Multiply a vector with the transposed matrix. y = At * x
223  void MultTranspose(const Vector &x, Vector &y) const;
224 
225  /// y += At * x (default) or y += a * At * x
226  void AddMultTranspose(const Vector &x, Vector &y,
227  const double a = 1.0) const;
228 
229  /** @brief Build and store internally the transpose of this matrix which will
230  be used in the methods AddMultTranspose() and MultTranspose(). */
231  /** If this method has been called, the internal transpose matrix will be
232  used to perform the action of the transpose matrix in AddMultTranspose(),
233  and MultTranspose().
234 
235  Warning: any changes in this matrix will invalidate the internal
236  transpose. To rebuild the transpose, call ResetTranspose() followed by a
237  call to this method. If the internal transpose is already built, this
238  method has no effect.
239 
240  When any non-default backend is enabled, i.e. Device::IsEnabled() is
241  true, the methods AddMultTranspose(), and MultTranspose(), require the
242  internal transpose to be built. If that is not the case (i.e. the
243  internal transpose is not built), these methods will raise an error with
244  an appropriate message pointing to this method. When using the default
245  backend, calling this method is optional.
246 
247  This method can only be used when the sparse matrix is finalized. */
248  void BuildTranspose() const;
249 
250  /** Reset (destroy) the internal transpose matrix. See BuildTranspose() for
251  more details. */
252  void ResetTranspose() const;
253 
254  void PartMult(const Array<int> &rows, const Vector &x, Vector &y) const;
255  void PartAddMult(const Array<int> &rows, const Vector &x, Vector &y,
256  const double a=1.0) const;
257 
258  /// y = A * x, treating all entries as booleans (zero=false, nonzero=true).
259  /** The actual values stored in the data array, #A, are not used - this means
260  and that all entries in the sparsity pattern are considered to be true by
261  this method. */
262  void BooleanMult(const Array<int> &x, Array<int> &y) const;
263 
264  /// y = At * x, treating all entries as booleans (zero=false, nonzero=true).
265  /** The actual values stored in the data array, #A, are not used - this means
266  and that all entries in the sparsity pattern are considered to be true by
267  this method. */
268  void BooleanMultTranspose(const Array<int> &x, Array<int> &y) const;
269 
270  /// Compute y^t A x
271  double InnerProduct(const Vector &x, const Vector &y) const;
272 
273  /// For all i compute \f$ x_i = \sum_j A_{ij} \f$
274  void GetRowSums(Vector &x) const;
275  /// For i = irow compute \f$ x_i = \sum_j | A_{i, j} | \f$
276  double GetRowNorml1(int irow) const;
277 
278  /// This virtual method is not supported: it always returns NULL.
279  virtual MatrixInverse *Inverse() const;
280 
281  /// Eliminates a column from the transpose matrix.
282  void EliminateRow(int row, const double sol, Vector &rhs);
283 
284  /// Eliminates a row from the matrix.
285  /*!
286  * - If @a dpolicy = #DIAG_ZERO, all the entries in the row will be set to 0.
287  * - If @a dpolicy = #DIAG_ONE (matrix must be square), the diagonal entry
288  * will be set equal to 1 and all other entries in the row to 0.
289  * - The policy #DIAG_KEEP is not supported.
290  */
291  void EliminateRow(int row, DiagonalPolicy dpolicy = DIAG_ZERO);
292 
293  /// Eliminates the column @a col from the matrix.
294  /** - If @a dpolicy = #DIAG_ZERO, all entries in the column will be set to 0.
295  - If @a dpolicy = #DIAG_ONE (matrix must be square), the diagonal entry
296  will be set equal to 1 and all other entries in the column to 0.
297  - The policy #DIAG_KEEP is not supported. */
298  void EliminateCol(int col, DiagonalPolicy dpolicy = DIAG_ZERO);
299 
300  /// Eliminate all columns i for which @a cols[i] != 0.
301  /** Elimination of a column means that all entries in the column are set to
302  zero. In addition, if the pointers @a x and @a b are not NULL, the
303  eliminated matrix entries are multiplied by the corresponding solution
304  value in @a *x and subtracted from the r.h.s. vector, @a *b. */
305  void EliminateCols(const Array<int> &cols, const Vector *x = NULL,
306  Vector *b = NULL);
307 
308  /// Eliminate row @a rc and column @a rc and modify the @a rhs using @a sol.
309  /** Eliminates the column @a rc to the @a rhs, deletes the row @a rc and
310  replaces the element (rc,rc) with 1.0; assumes that element (i,rc)
311  is assembled if and only if the element (rc,i) is assembled.
312  By default, elements (rc,rc) are set to 1.0, although this behavior
313  can be adjusted by changing the @a dpolicy parameter. */
314  void EliminateRowCol(int rc, const double sol, Vector &rhs,
315  DiagonalPolicy dpolicy = DIAG_ONE);
316 
317  /** @brief Similar to
318  EliminateRowCol(int, const double, Vector &, DiagonalPolicy), but
319  multiple values for eliminated unknowns are accepted, and accordingly
320  multiple right-hand-sides are used. */
321  void EliminateRowColMultipleRHS(int rc, const Vector &sol,
322  DenseMatrix &rhs,
323  DiagonalPolicy dpolicy = DIAG_ONE);
324 
325  /// Perform elimination and set the diagonal entry to the given value
326  void EliminateRowColDiag(int rc, double value);
327 
328  /// Eliminate row @a rc and column @a rc.
329  void EliminateRowCol(int rc, DiagonalPolicy dpolicy = DIAG_ONE);
330 
331  /** @brief Similar to EliminateRowCol(int, DiagonalPolicy) + save the
332  eliminated entries into @a Ae so that (*this) + Ae is equal to the
333  original matrix */
334  void EliminateRowCol(int rc, SparseMatrix &Ae,
335  DiagonalPolicy dpolicy = DIAG_ONE);
336 
337  /// If a row contains only one diag entry of zero, set it to 1.
338  void SetDiagIdentity();
339  /// If a row contains only zeros, set its diagonal to 1.
340  virtual void EliminateZeroRows(const double threshold = 1e-12);
341 
342  /// Gauss-Seidel forward and backward iterations over a vector x.
343  void Gauss_Seidel_forw(const Vector &x, Vector &y) const;
344  void Gauss_Seidel_back(const Vector &x, Vector &y) const;
345 
346  /// Determine appropriate scaling for Jacobi iteration
347  double GetJacobiScaling() const;
348  /** One scaled Jacobi iteration for the system A x = b.
349  x1 = x0 + sc D^{-1} (b - A x0) where D is the diag of A. */
350  void Jacobi(const Vector &b, const Vector &x0, Vector &x1, double sc) const;
351 
352  void DiagScale(const Vector &b, Vector &x, double sc = 1.0) const;
353 
354  /** x1 = x0 + sc D^{-1} (b - A x0) where \f$ D_{ii} = \sum_j |A_{ij}| \f$. */
355  void Jacobi2(const Vector &b, const Vector &x0, Vector &x1,
356  double sc = 1.0) const;
357 
358  /** x1 = x0 + sc D^{-1} (b - A x0) where \f$ D_{ii} = \sum_j A_{ij} \f$. */
359  void Jacobi3(const Vector &b, const Vector &x0, Vector &x1,
360  double sc = 1.0) const;
361 
362  /** @brief Finalize the matrix initialization, switching the storage format
363  from LIL to CSR. */
364  /** This method should be called once, after the matrix has been initialized.
365  Internally, this method converts the matrix from row-wise linked list
366  (LIL) format into CSR (compressed sparse row) format. */
367  virtual void Finalize(int skip_zeros = 1) { Finalize(skip_zeros, false); }
368 
369  /// A slightly more general version of the Finalize(int) method.
370  void Finalize(int skip_zeros, bool fix_empty_rows);
371 
372  bool Finalized() const { return !A.Empty(); }
373  bool areColumnsSorted() const { return isSorted; }
374 
375  /** @brief Remove entries smaller in absolute value than a given tolerance
376  @a tol. If @a fix_empty_rows is true, a zero value is inserted in the
377  diagonal entry (for square matrices only) */
378  void Threshold(double tol, bool fix_empty_rows = false);
379 
380  /** Split the matrix into M x N blocks of sparse matrices in CSR format.
381  The 'blocks' array is M x N (i.e. M and N are determined by its
382  dimensions) and its entries are overwritten by the new blocks. */
383  void GetBlocks(Array2D<SparseMatrix *> &blocks) const;
384 
385  void GetSubMatrix(const Array<int> &rows, const Array<int> &cols,
386  DenseMatrix &subm) const;
387 
388  /** @brief Initialize the SparseMatrix for fast access to the entries of the
389  given @a row which becomes the "current row". */
390  /** Fast access to the entries of the "current row" can be performed using
391  the methods: SearchRow(const int), _Add_(const int, const double),
392  _Set_(const int, const double), and _Get_(const int). */
393  inline void SetColPtr(const int row) const;
394  /** @brief Reset the "current row" set by calling SetColPtr(). This method
395  must be called between any two calls to SetColPtr(). */
396  inline void ClearColPtr() const;
397  /// Perform a fast search for an entry in the "current row". See SetColPtr().
398  /** If the matrix is not finalized and the entry is not found in the
399  SparseMatrix, it will be added to the sparsity pattern initialized with
400  zero. If the matrix is finalized and the entry is not found, an error
401  will be generated. */
402  inline double &SearchRow(const int col);
403  /// Add a value to an entry in the "current row". See SetColPtr().
404  inline void _Add_(const int col, const double a)
405  { SearchRow(col) += a; }
406  /// Set an entry in the "current row". See SetColPtr().
407  inline void _Set_(const int col, const double a)
408  { SearchRow(col) = a; }
409  /// Read the value of an entry in the "current row". See SetColPtr().
410  inline double _Get_(const int col) const;
411 
412  inline double &SearchRow(const int row, const int col);
413  inline void _Add_(const int row, const int col, const double a)
414  { SearchRow(row, col) += a; }
415  inline void _Set_(const int row, const int col, const double a)
416  { SearchRow(row, col) = a; }
417 
418  void Set(const int i, const int j, const double a);
419  void Add(const int i, const int j, const double a);
420 
421  void SetSubMatrix(const Array<int> &rows, const Array<int> &cols,
422  const DenseMatrix &subm, int skip_zeros = 1);
423 
424  void SetSubMatrixTranspose(const Array<int> &rows, const Array<int> &cols,
425  const DenseMatrix &subm, int skip_zeros = 1);
426 
427  void AddSubMatrix(const Array<int> &rows, const Array<int> &cols,
428  const DenseMatrix &subm, int skip_zeros = 1);
429 
430  bool RowIsEmpty(const int row) const;
431 
432  /// Extract all column indices and values from a given row.
433  /** If the matrix is finalized (i.e. in CSR format), @a cols and @a srow will
434  simply be references to the specific portion of the #J and #A arrays.
435  As required by the AbstractSparseMatrix interface this method returns:
436  - 0, if @a cols and @a srow are copies of the values in the matrix, i.e.
437  when the matrix is open.
438  - 1, if @a cols and @a srow are views of the values in the matrix, i.e.
439  when the matrix is finalized.
440  @warning This method breaks the const-ness when the matrix is finalized
441  because it gives write access to the #J and #A arrays. */
442  virtual int GetRow(const int row, Array<int> &cols, Vector &srow) const;
443 
444  void SetRow(const int row, const Array<int> &cols, const Vector &srow);
445  void AddRow(const int row, const Array<int> &cols, const Vector &srow);
446 
447  void ScaleRow(const int row, const double scale);
448  /// this = diag(sl) * this;
449  void ScaleRows(const Vector & sl);
450  /// this = this * diag(sr);
451  void ScaleColumns(const Vector & sr);
452 
453  /** @brief Add the sparse matrix 'B' to '*this'. This operation will cause an
454  error if '*this' is finalized and 'B' has larger sparsity pattern. */
456 
457  /** @brief Add the sparse matrix 'B' scaled by the scalar 'a' into '*this'.
458  Only entries in the sparsity pattern of '*this' are added. */
459  void Add(const double a, const SparseMatrix &B);
460 
461  SparseMatrix &operator=(double a);
462 
463  SparseMatrix &operator*=(double a);
464 
465  /// Prints matrix to stream out.
466  void Print(std::ostream &out = mfem::out, int width_ = 4) const;
467 
468  /// Prints matrix in matlab format.
469  void PrintMatlab(std::ostream &out = mfem::out) const;
470 
471  /// Prints matrix in Matrix Market sparse format.
472  void PrintMM(std::ostream &out = mfem::out) const;
473 
474  /// Prints matrix to stream out in hypre_CSRMatrix format.
475  void PrintCSR(std::ostream &out) const;
476 
477  /// Prints a sparse matrix to stream out in CSR format.
478  void PrintCSR2(std::ostream &out) const;
479 
480  /// Print various sparse matrix staticstics.
481  void PrintInfo(std::ostream &out) const;
482 
483  /// Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix
484  double IsSymmetric() const;
485 
486  /// (*this) = 1/2 ((*this) + (*this)^t)
487  void Symmetrize();
488 
489  /// Returns the number of the nonzero elements in the matrix
490  virtual int NumNonZeroElems() const;
491 
492  double MaxNorm() const;
493 
494  /// Count the number of entries with |a_ij| <= tol.
495  int CountSmallElems(double tol) const;
496 
497  /// Count the number of entries that are NOT finite, i.e. Inf or Nan.
498  int CheckFinite() const;
499 
500  /// Set the graph ownership flag (I and J arrays).
501  void SetGraphOwner(bool ownij)
502  { I.SetHostPtrOwner(ownij); J.SetHostPtrOwner(ownij); }
503 
504  /// Set the data ownership flag (A array).
505  void SetDataOwner(bool owna) { A.SetHostPtrOwner(owna); }
506 
507  /// Get the graph ownership flag (I and J arrays).
508  bool OwnsGraph() const { return I.OwnsHostPtr() && J.OwnsHostPtr(); }
509 
510  /// Get the data ownership flag (A array).
511  bool OwnsData() const { return A.OwnsHostPtr(); }
512 
513  /// Lose the ownership of the graph (I, J) and data (A) arrays.
514  void LoseData() { SetGraphOwner(false); SetDataOwner(false); }
515 
516  void Swap(SparseMatrix &other);
517 
518  /// Destroys sparse matrix.
519  virtual ~SparseMatrix() { Destroy(); }
520 
521  Type GetType() const { return MFEM_SPARSEMAT; }
522 };
523 
524 /// Applies f() to each element of the matrix (after it is finalized).
525 void SparseMatrixFunction(SparseMatrix &S, double (*f)(double));
526 
527 
528 /// Transpose of a sparse matrix. A must be finalized.
529 SparseMatrix *Transpose(const SparseMatrix &A);
530 /// Transpose of a sparse matrix. A does not need to be a CSR matrix.
531 SparseMatrix *TransposeAbstractSparseMatrix (const AbstractSparseMatrix &A,
532  int useActualWidth);
533 
534 /// Matrix product A.B.
535 /** If @a OAB is not NULL, we assume it has the structure of A.B and store the
536  result in @a OAB. If @a OAB is NULL, we create a new SparseMatrix to store
537  the result and return a pointer to it.
538 
539  All matrices must be finalized. */
540 SparseMatrix *Mult(const SparseMatrix &A, const SparseMatrix &B,
541  SparseMatrix *OAB = NULL);
542 
543 /// C = A^T B
544 SparseMatrix *TransposeMult(const SparseMatrix &A, const SparseMatrix &B);
545 
546 /// Matrix product of sparse matrices. A and B do not need to be CSR matrices
547 SparseMatrix *MultAbstractSparseMatrix (const AbstractSparseMatrix &A,
548  const AbstractSparseMatrix &B);
549 
550 /// Matrix product A.B
551 DenseMatrix *Mult(const SparseMatrix &A, DenseMatrix &B);
552 
553 /// RAP matrix product (with R=P^T)
554 DenseMatrix *RAP(const SparseMatrix &A, DenseMatrix &P);
555 
556 /// RAP matrix product (with R=P^T)
557 DenseMatrix *RAP(DenseMatrix &A, const SparseMatrix &P);
558 
559 /** RAP matrix product (with P=R^T). ORAP is like OAB above.
560  All matrices must be finalized. */
561 SparseMatrix *RAP(const SparseMatrix &A, const SparseMatrix &R,
562  SparseMatrix *ORAP = NULL);
563 
564 /// General RAP with given R^T, A and P
565 SparseMatrix *RAP(const SparseMatrix &Rt, const SparseMatrix &A,
566  const SparseMatrix &P);
567 
568 /// Matrix multiplication A^t D A. All matrices must be finalized.
569 SparseMatrix *Mult_AtDA(const SparseMatrix &A, const Vector &D,
570  SparseMatrix *OAtDA = NULL);
571 
572 
573 /// Matrix addition result = A + B.
574 SparseMatrix * Add(const SparseMatrix & A, const SparseMatrix & B);
575 /// Matrix addition result = a*A + b*B
576 SparseMatrix * Add(double a, const SparseMatrix & A, double b,
577  const SparseMatrix & B);
578 /// Matrix addition result = sum_i A_i
579 SparseMatrix * Add(Array<SparseMatrix *> & Ai);
580 
581 /// B += alpha * A
582 void Add(const SparseMatrix &A, double alpha, DenseMatrix &B);
583 
584 /// Produces a block matrix with blocks A_{ij}*B
585 DenseMatrix *OuterProduct(const DenseMatrix &A, const DenseMatrix &B);
586 
587 /// Produces a block matrix with blocks A_{ij}*B
588 SparseMatrix *OuterProduct(const DenseMatrix &A, const SparseMatrix &B);
589 
590 /// Produces a block matrix with blocks A_{ij}*B
591 SparseMatrix *OuterProduct(const SparseMatrix &A, const DenseMatrix &B);
592 
593 /// Produces a block matrix with blocks A_{ij}*B
594 SparseMatrix *OuterProduct(const SparseMatrix &A, const SparseMatrix &B);
595 
596 
597 // Inline methods
598 
599 inline void SparseMatrix::SetColPtr(const int row) const
600 {
601  if (Rows)
602  {
603  if (ColPtrNode == NULL)
604  {
605  ColPtrNode = new RowNode *[width];
606  for (int i = 0; i < width; i++)
607  {
608  ColPtrNode[i] = NULL;
609  }
610  }
611  for (RowNode *node_p = Rows[row]; node_p != NULL; node_p = node_p->Prev)
612  {
613  ColPtrNode[node_p->Column] = node_p;
614  }
615  }
616  else
617  {
618  if (ColPtrJ == NULL)
619  {
620  ColPtrJ = new int[width];
621  for (int i = 0; i < width; i++)
622  {
623  ColPtrJ[i] = -1;
624  }
625  }
626  for (int j = I[row], end = I[row+1]; j < end; j++)
627  {
628  ColPtrJ[J[j]] = j;
629  }
630  }
631  current_row = row;
632 }
633 
634 inline void SparseMatrix::ClearColPtr() const
635 {
636  if (Rows)
637  {
638  for (RowNode *node_p = Rows[current_row]; node_p != NULL;
639  node_p = node_p->Prev)
640  {
641  ColPtrNode[node_p->Column] = NULL;
642  }
643  }
644  else
645  {
646  for (int j = I[current_row], end = I[current_row+1]; j < end; j++)
647  {
648  ColPtrJ[J[j]] = -1;
649  }
650  }
651 }
652 
653 inline double &SparseMatrix::SearchRow(const int col)
654 {
655  if (Rows)
656  {
657  RowNode *node_p = ColPtrNode[col];
658  if (node_p == NULL)
659  {
660 #ifdef MFEM_USE_MEMALLOC
661  node_p = NodesMem->Alloc();
662 #else
663  node_p = new RowNode;
664 #endif
665  node_p->Prev = Rows[current_row];
666  node_p->Column = col;
667  node_p->Value = 0.0;
668  Rows[current_row] = ColPtrNode[col] = node_p;
669  }
670  return node_p->Value;
671  }
672  else
673  {
674  const int j = ColPtrJ[col];
675  MFEM_VERIFY(j != -1, "Entry for column " << col << " is not allocated.");
676  return A[j];
677  }
678 }
679 
680 inline double SparseMatrix::_Get_(const int col) const
681 {
682  if (Rows)
683  {
684  RowNode *node_p = ColPtrNode[col];
685  return (node_p == NULL) ? 0.0 : node_p->Value;
686  }
687  else
688  {
689  const int j = ColPtrJ[col];
690  return (j == -1) ? 0.0 : A[j];
691  }
692 }
693 
694 inline double &SparseMatrix::SearchRow(const int row, const int col)
695 {
696  if (Rows)
697  {
698  RowNode *node_p;
699 
700  for (node_p = Rows[row]; 1; node_p = node_p->Prev)
701  {
702  if (node_p == NULL)
703  {
704 #ifdef MFEM_USE_MEMALLOC
705  node_p = NodesMem->Alloc();
706 #else
707  node_p = new RowNode;
708 #endif
709  node_p->Prev = Rows[row];
710  node_p->Column = col;
711  node_p->Value = 0.0;
712  Rows[row] = node_p;
713  break;
714  }
715  else if (node_p->Column == col)
716  {
717  break;
718  }
719  }
720  return node_p->Value;
721  }
722  else
723  {
724  int *Ip = I+row, *Jp = J;
725  for (int k = Ip[0], end = Ip[1]; k < end; k++)
726  {
727  if (Jp[k] == col)
728  {
729  return A[k];
730  }
731  }
732  MFEM_ABORT("Could not find entry for row = " << row << ", col = " << col);
733  }
734  return A[0];
735 }
736 
737 /// Specialization of the template function Swap<> for class SparseMatrix
738 template<> inline void Swap<SparseMatrix>(SparseMatrix &a, SparseMatrix &b)
739 {
740  a.Swap(b);
741 }
742 
743 } // namespace mfem
744 
745 #endif
RowNode ** ColPtrNode
Definition: sparsemat.hpp:67
Memory< int > I
Array with size (height+1) containing the row offsets.
Definition: sparsemat.hpp:52
Elem * Alloc()
Definition: mem_alloc.hpp:146
void _Add_(const int row, const int col, const double a)
Definition: sparsemat.hpp:413
double GetJacobiScaling() const
Determine appropriate scaling for Jacobi iteration.
Definition: sparsemat.cpp:2043
int CheckFinite() const
Count the number of entries that are NOT finite, i.e. Inf or Nan.
Definition: sparsemat.cpp:1234
static MemoryClass GetMemoryClass()
Get the current Device MemoryClass. This is the MemoryClass used by most MFEM device kernels to acces...
Definition: device.hpp:214
SparseMatrix * TransposeAbstractSparseMatrix(const AbstractSparseMatrix &A, int useActualWidth)
Transpose of a sparse matrix. A does not need to be a CSR matrix.
Definition: sparsemat.cpp:2997
int RowSize(const int i) const
Returns the number of elements in row i.
Definition: sparsemat.cpp:252
virtual int NumNonZeroElems() const
Returns the number of the nonzero elements in the matrix.
Definition: sparsemat.cpp:1159
void Add(const int i, const int j, const double a)
Definition: sparsemat.cpp:2244
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:404
void BuildTranspose() const
Build and store internally the transpose of this matrix which will be used in the methods AddMultTran...
Definition: sparsemat.cpp:668
void Clear()
Clear the contents of the SparseMatrix.
Definition: sparsemat.hpp:135
void EliminateCol(int col, DiagonalPolicy dpolicy=DIAG_ZERO)
Eliminates the column col from the matrix.
Definition: sparsemat.cpp:1311
SparseMatrix * At
Transpose of A. Owned. Used to perform MultTranspose() on devices.
Definition: sparsemat.hpp:70
Type GetType() const
Definition: sparsemat.hpp:521
void _Set_(const int row, const int col, const double a)
Definition: sparsemat.hpp:415
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1563
void DiagScale(const Vector &b, Vector &x, double sc=1.0) const
Definition: sparsemat.cpp:2107
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
Definition: sparsemat.hpp:367
void MakeRef(const SparseMatrix &master)
Clear the contents of the SparseMatrix and make it a reference to master.
Definition: sparsemat.cpp:223
double & SearchRow(const int col)
Perform a fast search for an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:653
DenseMatrix * ToDenseMatrix() const
Produces a DenseMatrix from a SparseMatrix.
Definition: sparsemat.cpp:515
bool Empty() const
Check if the SparseMatrix is empty.
Definition: sparsemat.hpp:138
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Definition: sparsemat.cpp:299
Set the diagonal value to zero.
Definition: matrix.hpp:34
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:726
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:146
Abstract data type for sparse matrices.
Definition: matrix.hpp:77
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:2586
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:552
double MaxNorm() const
Definition: sparsemat.cpp:1181
virtual ~SparseMatrix()
Destroys sparse matrix.
Definition: sparsemat.hpp:519
int * GetI()
Return the array I.
Definition: sparsemat.hpp:141
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:623
Abstract data type for matrix inverse.
Definition: matrix.hpp:66
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:1486
void PrintInfo(std::ostream &out) const
Print various sparse matrix staticstics.
Definition: sparsemat.cpp:2822
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
Definition: sparsemat.cpp:418
void LoseData()
Lose the ownership of the graph (I, J) and data (A) arrays.
Definition: sparsemat.hpp:514
void AddRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:2467
double * GetRowEntries(const int row)
Return a pointer to the entries in a row.
Definition: sparsemat.cpp:313
void SortColumnIndices()
Sort the column indices corresponding to each row.
Definition: sparsemat.cpp:365
void PrintMM(std::ostream &out=mfem::out) const
Prints matrix in Matrix Market sparse format.
Definition: sparsemat.cpp:2751
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const
Definition: sparsemat.cpp:2331
void _Set_(const int col, const double a)
Set an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:407
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:756
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:3365
void PartAddMult(const Array< int > &rows, const Vector &x, Vector &y, const double a=1.0) const
Definition: sparsemat.cpp:708
bool RowIsEmpty(const int row) const
Definition: sparsemat.cpp:2359
void GetBlocks(Array2D< SparseMatrix * > &blocks) const
Definition: sparsemat.cpp:1008
void ScaleRow(const int row, const double scale)
Definition: sparsemat.cpp:2499
RowNodeAlloc * NodesMem
Definition: sparsemat.hpp:74
double * GetData()
Return the element data, i.e. the array A.
Definition: sparsemat.hpp:151
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
Definition: sparsemat.cpp:1140
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:394
void SetGraphOwner(bool ownij)
Set the graph ownership flag (I and J arrays).
Definition: sparsemat.hpp:501
ID for class SparseMatrix.
Definition: operator.hpp:138
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:3030
bool OwnsData() const
Get the data ownership flag (A array).
Definition: sparsemat.hpp:511
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:634
void ScaleColumns(const Vector &sr)
this = this * diag(sr);
Definition: sparsemat.cpp:2558
void PrintCSR2(std::ostream &out) const
Prints a sparse matrix to stream out in CSR format.
Definition: sparsemat.cpp:2797
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:2380
bool isSorted
Are the columns sorted already.
Definition: sparsemat.hpp:78
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:58
SparseMatrix()
Create an empty SparseMatrix.
Definition: sparsemat.hpp:85
Data type sparse matrix.
Definition: sparsemat.hpp:40
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:36
bool OwnsGraph() const
Get the graph ownership flag (I and J arrays).
Definition: sparsemat.hpp:508
void SetWidth(int width_=-1)
Change the width of a SparseMatrix.
Definition: sparsemat.cpp:327
bool areColumnsSorted() const
Definition: sparsemat.hpp:373
double IsSymmetric() const
Returns max_{i,j} |(i,j)-(j,i)| for a finalized matrix.
Definition: sparsemat.cpp:1098
void SetDiagIdentity()
If a row contains only one diag entry of zero, set it to 1.
Definition: sparsemat.cpp:1846
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
Definition: table.cpp:414
bool OwnsHostPtr() const
Return true if the host pointer is owned. Ownership indicates whether the pointer will be deleted by ...
void ScaleRows(const Vector &sl)
this = diag(sl) * this;
Definition: sparsemat.cpp:2527
void Swap< SparseMatrix >(SparseMatrix &a, SparseMatrix &b)
Specialization of the template function Swap&lt;&gt; for class SparseMatrix.
Definition: sparsemat.hpp:738
DenseMatrix * OuterProduct(const DenseMatrix &A, const DenseMatrix &B)
Produces a block matrix with blocks A_{ij}*B.
Definition: sparsemat.cpp:3509
void PrintMatlab(std::ostream &out=mfem::out) const
Prints matrix in matlab format.
Definition: sparsemat.cpp:2731
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2187
Set the diagonal value to one.
Definition: matrix.hpp:35
Dynamic 2D array using row-major layout.
Definition: array.hpp:304
void Jacobi(const Vector &b, const Vector &x0, Vector &x1, double sc) const
Definition: sparsemat.cpp:2076
void SetSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
Definition: sparsemat.cpp:2263
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
Definition: sparsemat.hpp:213
bool Finalized() const
Definition: sparsemat.hpp:372
const double * GetData() const
Return the element data, i.e. the array A, const version.
Definition: sparsemat.hpp:153
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:135
void SparseMatrixFunction(SparseMatrix &S, double(*f)(double))
Applies f() to each element of the matrix (after it is finalized).
Definition: sparsemat.cpp:2930
virtual MatrixInverse * Inverse() const
This virtual method is not supported: it always returns NULL.
Definition: sparsemat.cpp:1258
void EliminateRowColDiag(int rc, double value)
Perform elimination and set the diagonal entry to the given value.
Definition: sparsemat.cpp:1689
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
Definition: sparsemat.cpp:545
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:599
const int * GetJ() const
Return the array J, const version.
Definition: sparsemat.hpp:148
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:2296
void Set(const int i, const int j, const double a)
Definition: sparsemat.cpp:2225
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:1390
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:867
int Size() const
For backward compatibility, define Size() to be synonym of Height().
Definition: sparsemat.hpp:132
void Gauss_Seidel_back(const Vector &x, Vector &y) const
Definition: sparsemat.cpp:1960
void Gauss_Seidel_forw(const Vector &x, Vector &y) const
Gauss-Seidel forward and backward iterations over a vector x.
Definition: sparsemat.cpp:1876
RowNode ** Rows
Array of linked lists, one for every row. This array represents the linked list (LIL) storage format...
Definition: sparsemat.hpp:63
void EliminateRow(int row, const double sol, Vector &rhs)
Eliminates a column from the transpose matrix.
Definition: sparsemat.cpp:1263
int ActualWidth() const
Returns the actual Width of the matrix.
Definition: sparsemat.cpp:2904
void GetDiag(Vector &d) const
Returns the Diagonal of A.
Definition: sparsemat.cpp:486
MemAlloc< RowNode, 1024 > RowNodeAlloc
Definition: sparsemat.hpp:73
Host memory; using new[] and delete[].
double GetRowNorml1(int irow) const
For i = irow compute .
Definition: sparsemat.cpp:843
void ResetTranspose() const
Definition: sparsemat.cpp:676
SparseMatrix * TransposeMult(const SparseMatrix &A, const SparseMatrix &B)
C = A^T B.
Definition: sparsemat.cpp:3202
virtual void EliminateZeroRows(const double threshold=1e-12)
If a row contains only zeros, set its diagonal to 1.
Definition: sparsemat.cpp:1857
bool Empty() const
Return true if the Memory object is empty, see Reset().
double & operator()(int i, int j)
Returns reference to A[i][j].
Definition: sparsemat.cpp:428
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:616
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:1351
SparseMatrix & operator=(const SparseMatrix &rhs)
Assignment operator: deep copy.
Definition: sparsemat.cpp:213
SparseMatrix & operator*=(double a)
Definition: sparsemat.cpp:2663
void SetDataOwner(bool owna)
Set the data ownership flag (A array).
Definition: sparsemat.hpp:505
void SetRow(const int row, const Array< int > &cols, const Vector &srow)
Definition: sparsemat.cpp:2420
const double alpha
Definition: ex15.cpp:339
Vector data type.
Definition: vector.hpp:48
void Print(std::ostream &out=mfem::out, int width_=4) const
Prints matrix to stream out.
Definition: sparsemat.cpp:2687
const int * GetI() const
Return the array I, const version.
Definition: sparsemat.hpp:143
void Jacobi2(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
Definition: sparsemat.cpp:2139
int MaxRowSize() const
Returns the maximum number of elements among all rows.
Definition: sparsemat.cpp:275
void Jacobi3(const Vector &b, const Vector &x0, Vector &x1, double sc=1.0) const
Definition: sparsemat.cpp:2163
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:3210
double _Get_(const int col) const
Read the value of an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:680
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:64
Memory< int > J
Array with size I[height], containing the column indices for all matrix entries, as indexed by the I ...
Definition: sparsemat.hpp:55
void Swap(SparseMatrix &other)
Definition: sparsemat.cpp:3617
void PrintCSR(std::ostream &out) const
Prints matrix to stream out in hypre_CSRMatrix format.
Definition: sparsemat.cpp:2773
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: sparsemat.cpp:779
MemoryClass
Memory classes identify subsets of memory types.
Definition: mem_manager.hpp:40
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:25
int CountSmallElems(double tol) const
Count the number of entries with |a_ij| &lt;= tol.
Definition: sparsemat.cpp:1206
void PartMult(const Array< int > &rows, const Vector &x, Vector &y) const
Definition: sparsemat.cpp:682
void GetRowSums(Vector &x) const
For all i compute .
Definition: sparsemat.cpp:820