MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
hypre.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_HYPRE
13 #define MFEM_HYPRE
14 
15 #include "../config/config.hpp"
16 
17 #ifdef MFEM_USE_MPI
18 
19 #include <mpi.h>
20 
21 // Enable internal hypre timing routines
22 #define HYPRE_TIMING
23 
24 // hypre header files
25 #include "seq_mv.h"
26 #include "_hypre_parcsr_mv.h"
27 #include "_hypre_parcsr_ls.h"
28 #include "temp_multivector.h"
29 #include "../general/globals.hpp"
30 
31 #ifdef HYPRE_COMPLEX
32 #error "MFEM does not work with HYPRE's complex numbers support"
33 #endif
34 
35 #include "sparsemat.hpp"
36 #include "hypre_parcsr.hpp"
37 
38 namespace mfem
39 {
40 
41 class ParFiniteElementSpace;
42 class HypreParMatrix;
43 
44 namespace internal
45 {
46 
47 template <typename int_type>
48 inline int to_int(int_type i)
49 {
50  MFEM_ASSERT(int_type(int(i)) == i, "overflow converting int_type to int");
51  return int(i);
52 }
53 
54 // Specialization for to_int(int)
55 template <> inline int to_int(int i) { return i; }
56 
57 // Convert a HYPRE_Int to int
58 #ifdef HYPRE_BIGINT
59 template <>
60 inline int to_int(HYPRE_Int i)
61 {
62  MFEM_ASSERT(HYPRE_Int(int(i)) == i, "overflow converting HYPRE_Int to int");
63  return int(i);
64 }
65 #endif
66 
67 }
68 
69 /// Wrapper for hypre's parallel vector class
70 class HypreParVector : public Vector
71 {
72 private:
73  int own_ParVector;
74 
75  /// The actual object
76  hypre_ParVector *x;
77 
78  friend class HypreParMatrix;
79 
80  // Set Vector::data and Vector::size from *x
81  inline void _SetDataAndSize_();
82 
83 public:
84  /** @brief Creates vector with given global size and parallel partitioning of
85  the rows/columns given by @a col. */
86  /** @anchor hypre_partitioning_descr
87  The partitioning is defined in one of two ways depending on the
88  configuration of HYPRE:
89  1. If HYPRE_AssumedPartitionCheck() returns true (the default),
90  then col is of length 2 and the local processor owns columns
91  [col[0],col[1]).
92  2. If HYPRE_AssumedPartitionCheck() returns false, then col is of
93  length (number of processors + 1) and processor P owns columns
94  [col[P],col[P+1]) i.e. each processor has a copy of the same col
95  array. */
96  HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col);
97  /** @brief Creates vector with given global size, partitioning of the
98  columns, and data. */
99  /** The data must be allocated and destroyed outside. If @a _data is NULL, a
100  dummy vector without a valid data array will be created. See @ref
101  hypre_partitioning_descr "here" for a description of the @a col array. */
102  HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, double *_data,
103  HYPRE_Int *col);
104  /// Creates vector compatible with y
105  HypreParVector(const HypreParVector &y);
106  /// Creates vector compatible with (i.e. in the domain of) A or A^T
107  explicit HypreParVector(const HypreParMatrix &A, int transpose = 0);
108  /// Creates vector wrapping y
109  explicit HypreParVector(HYPRE_ParVector y);
110  /// Create a true dof parallel vector on a given ParFiniteElementSpace
111  explicit HypreParVector(ParFiniteElementSpace *pfes);
112 
113  /// MPI communicator
114  MPI_Comm GetComm() { return x->comm; }
115 
116  /// Returns the parallel row/column partitioning
117  /** See @ref hypre_partitioning_descr "here" for a description of the
118  partitioning array. */
119  inline HYPRE_Int *Partitioning() { return x->partitioning; }
120 
121  /// Returns the global number of rows
122  inline HYPRE_Int GlobalSize() { return x->global_size; }
123 
124  /// Typecasting to hypre's hypre_ParVector*
125  operator hypre_ParVector*() const { return x; }
126 #ifndef HYPRE_PAR_VECTOR_STRUCT
127  /// Typecasting to hypre's HYPRE_ParVector, a.k.a. void *
128  operator HYPRE_ParVector() const { return (HYPRE_ParVector) x; }
129 #endif
130  /// Changes the ownership of the the vector
131  hypre_ParVector *StealParVector() { own_ParVector = 0; return x; }
132 
133  /// Sets ownership of the internal hypre_ParVector
134  void SetOwnership(int own) { own_ParVector = own; }
135 
136  /// Gets ownership of the internal hypre_ParVector
137  int GetOwnership() const { return own_ParVector; }
138 
139  /// Returns the global vector in each processor
140  Vector* GlobalVector() const;
141 
142  /// Set constant values
143  HypreParVector& operator= (double d);
144  /// Define '=' for hypre vectors.
146 
147  /// Sets the data of the Vector and the hypre_ParVector to @a _data.
148  /** Must be used only for HypreParVector%s that do not own the data,
149  e.g. created with the constructor:
150  HypreParVector(MPI_Comm, HYPRE_Int, double *, HYPRE_Int *). */
151  void SetData(double *_data);
152 
153  /// Set random values
154  HYPRE_Int Randomize(HYPRE_Int seed);
155 
156  /// Prints the locally owned rows in parallel
157  void Print(const char *fname) const;
158 
159  /// Calls hypre's destroy function
160  ~HypreParVector();
161 
162 #ifdef MFEM_USE_SUNDIALS
163  /// (DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_PARALLEL.
164  /** @deprecated The returned N_Vector must be destroyed by the caller. */
165  MFEM_DEPRECATED virtual N_Vector ToNVector();
166  using Vector::ToNVector;
167 #endif
168 };
169 
170 /// Returns the inner product of x and y
171 double InnerProduct(HypreParVector &x, HypreParVector &y);
172 double InnerProduct(HypreParVector *x, HypreParVector *y);
173 
174 
175 /** @brief Compute the l_p norm of the Vector which is split without overlap
176  across the given communicator. */
177 double ParNormlp(const Vector &vec, double p, MPI_Comm comm);
178 
179 
180 /// Wrapper for hypre's ParCSR matrix class
181 class HypreParMatrix : public Operator
182 {
183 private:
184  /// The actual object
185  hypre_ParCSRMatrix *A;
186 
187  /// Auxiliary vectors for typecasting
188  mutable HypreParVector *X, *Y;
189 
190  // Flags indicating ownership of A->diag->{i,j,data}, A->offd->{i,j,data},
191  // and A->col_map_offd.
192  // The possible values for diagOwner are:
193  // -1: no special treatment of A->diag (default)
194  // 0: prevent hypre from destroying A->diag->{i,j,data}
195  // 1: same as 0, plus take ownership of A->diag->{i,j}
196  // 2: same as 0, plus take ownership of A->diag->data
197  // 3: same as 0, plus take ownership of A->diag->{i,j,data}
198  // The same values and rules apply to offdOwner and A->offd.
199  // The possible values for colMapOwner are:
200  // -1: no special treatment of A->col_map_offd (default)
201  // 0: prevent hypre from destroying A->col_map_offd
202  // 1: same as 0, plus take ownership of A->col_map_offd
203  // All owned arrays are destroyed with 'delete []'.
204  signed char diagOwner, offdOwner, colMapOwner;
205 
206  // Does the object own the pointer A?
207  signed char ParCSROwner;
208 
209  // Initialize with defaults. Does not initialize inherited members.
210  void Init();
211 
212  // Delete all owned data. Does not perform re-initialization with defaults.
213  void Destroy();
214 
215  // Copy (shallow/deep, based on HYPRE_BIGINT) the I and J arrays from csr to
216  // hypre_csr. Shallow copy the data. Return the appropriate ownership flag.
217  static char CopyCSR(SparseMatrix *csr, hypre_CSRMatrix *hypre_csr);
218  // Copy (shallow or deep, based on HYPRE_BIGINT) the I and J arrays from
219  // bool_csr to hypre_csr. Allocate the data array and set it to all ones.
220  // Return the appropriate ownership flag.
221  static char CopyBoolCSR(Table *bool_csr, hypre_CSRMatrix *hypre_csr);
222 
223  // Copy the j array of a hypre_CSRMatrix to the given J array, converting
224  // the indices from HYPRE_Int to int.
225  static void CopyCSR_J(hypre_CSRMatrix *hypre_csr, int *J);
226 
227 public:
228  /// An empty matrix to be used as a reference to an existing matrix
229  HypreParMatrix();
230 
231  /// Converts hypre's format to HypreParMatrix
232  /** If @a owner is false, ownership of @a a is not transferred */
233  explicit HypreParMatrix(hypre_ParCSRMatrix *a, bool owner = true)
234  {
235  Init();
236  A = a;
237  if (!owner) { ParCSROwner = 0; }
238  height = GetNumRows();
239  width = GetNumCols();
240  }
241 
242  /// Creates block-diagonal square parallel matrix.
243  /** Diagonal is given by @a diag which must be in CSR format (finalized). The
244  new HypreParMatrix does not take ownership of any of the input arrays.
245  See @ref hypre_partitioning_descr "here" for a description of the row
246  partitioning array @a row_starts.
247 
248  @warning The ordering of the columns in each row in @a *diag may be
249  changed by this constructor to ensure that the first entry in each row is
250  the diagonal one. This is expected by most hypre functions. */
251  HypreParMatrix(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *row_starts,
252  SparseMatrix *diag); // constructor with 4 arguments, v1
253 
254  /// Creates block-diagonal rectangular parallel matrix.
255  /** Diagonal is given by @a diag which must be in CSR format (finalized). The
256  new HypreParMatrix does not take ownership of any of the input arrays.
257  See @ref hypre_partitioning_descr "here" for a description of the
258  partitioning arrays @a row_starts and @a col_starts. */
259  HypreParMatrix(MPI_Comm comm, HYPRE_Int global_num_rows,
260  HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
261  HYPRE_Int *col_starts,
262  SparseMatrix *diag); // constructor with 6 arguments, v1
263 
264  /// Creates general (rectangular) parallel matrix.
265  /** The new HypreParMatrix does not take ownership of any of the input
266  arrays. See @ref hypre_partitioning_descr "here" for a description of the
267  partitioning arrays @a row_starts and @a col_starts. */
268  HypreParMatrix(MPI_Comm comm, HYPRE_Int global_num_rows,
269  HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
270  HYPRE_Int *col_starts, SparseMatrix *diag, SparseMatrix *offd,
271  HYPRE_Int *cmap); // constructor with 8 arguments
272 
273  /// Creates general (rectangular) parallel matrix.
274  /** The new HypreParMatrix takes ownership of all input arrays, except
275  @a col_starts and @a row_starts. See @ref hypre_partitioning_descr "here"
276  for a description of the partitioning arrays @a row_starts and @a
277  col_starts. */
278  HypreParMatrix(MPI_Comm comm,
279  HYPRE_Int global_num_rows, HYPRE_Int global_num_cols,
280  HYPRE_Int *row_starts, HYPRE_Int *col_starts,
281  HYPRE_Int *diag_i, HYPRE_Int *diag_j, double *diag_data,
282  HYPRE_Int *offd_i, HYPRE_Int *offd_j, double *offd_data,
283  HYPRE_Int offd_num_cols,
284  HYPRE_Int *offd_col_map); // constructor with 13 arguments
285 
286  /// Creates a parallel matrix from SparseMatrix on processor 0.
287  /** See @ref hypre_partitioning_descr "here" for a description of the
288  partitioning arrays @a row_starts and @a col_starts. */
289  HypreParMatrix(MPI_Comm comm, HYPRE_Int *row_starts, HYPRE_Int *col_starts,
290  SparseMatrix *a); // constructor with 4 arguments, v2
291 
292  /// Creates boolean block-diagonal rectangular parallel matrix.
293  /** The new HypreParMatrix does not take ownership of any of the input
294  arrays. See @ref hypre_partitioning_descr "here" for a description of the
295  partitioning arrays @a row_starts and @a col_starts. */
296  HypreParMatrix(MPI_Comm comm, HYPRE_Int global_num_rows,
297  HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
298  HYPRE_Int *col_starts,
299  Table *diag); // constructor with 6 arguments, v2
300 
301  /// Creates boolean rectangular parallel matrix.
302  /** The new HypreParMatrix takes ownership of the arrays @a i_diag,
303  @a j_diag, @a i_offd, @a j_offd, and @a cmap; does not take ownership of
304  the arrays @a row and @a col. See @ref hypre_partitioning_descr "here"
305  for a description of the partitioning arrays @a row and @a col. */
306  HypreParMatrix(MPI_Comm comm, int id, int np, HYPRE_Int *row, HYPRE_Int *col,
307  HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd,
308  HYPRE_Int *j_offd, HYPRE_Int *cmap,
309  HYPRE_Int cmap_size); // constructor with 11 arguments
310 
311  /** @brief Creates a general parallel matrix from a local CSR matrix on each
312  processor described by the @a I, @a J and @a data arrays. */
313  /** The local matrix should be of size (local) @a nrows by (global)
314  @a glob_ncols. The new parallel matrix contains copies of all input
315  arrays (so they can be deleted). See @ref hypre_partitioning_descr "here"
316  for a description of the partitioning arrays @a rows and @a cols. */
317  HypreParMatrix(MPI_Comm comm, int nrows, HYPRE_Int glob_nrows,
318  HYPRE_Int glob_ncols, int *I, HYPRE_Int *J,
319  double *data, HYPRE_Int *rows,
320  HYPRE_Int *cols); // constructor with 9 arguments
321 
322  /** @brief Copy constructor for a ParCSR matrix which creates a deep copy of
323  structure and data from @a P. */
324  HypreParMatrix(const HypreParMatrix &P);
325 
326  /// Make this HypreParMatrix a reference to 'master'
327  void MakeRef(const HypreParMatrix &master);
328 
329  /// MPI communicator
330  MPI_Comm GetComm() const { return A->comm; }
331 
332  /// Typecasting to hypre's hypre_ParCSRMatrix*
333  operator hypre_ParCSRMatrix*() const { return A; }
334 #ifndef HYPRE_PAR_CSR_MATRIX_STRUCT
335  /// Typecasting to hypre's HYPRE_ParCSRMatrix, a.k.a. void *
336  operator HYPRE_ParCSRMatrix() { return (HYPRE_ParCSRMatrix) A; }
337 #endif
338  /// Changes the ownership of the the matrix
339  hypre_ParCSRMatrix* StealData();
340 
341  /// Explicitly set the three ownership flags, see docs for diagOwner etc.
342  void SetOwnerFlags(signed char diag, signed char offd, signed char colmap)
343  { diagOwner = diag, offdOwner = offd, colMapOwner = colmap; }
344 
345  /// Get diag ownership flag
346  signed char OwnsDiag() const { return diagOwner; }
347  /// Get offd ownership flag
348  signed char OwnsOffd() const { return offdOwner; }
349  /// Get colmap ownership flag
350  signed char OwnsColMap() const { return colMapOwner; }
351 
352  /** If the HypreParMatrix does not own the row-starts array, make a copy of
353  it that the HypreParMatrix will own. If the col-starts array is the same
354  as the row-starts array, col-starts is also replaced. */
355  void CopyRowStarts();
356  /** If the HypreParMatrix does not own the col-starts array, make a copy of
357  it that the HypreParMatrix will own. If the row-starts array is the same
358  as the col-starts array, row-starts is also replaced. */
359  void CopyColStarts();
360 
361  /// Returns the global number of nonzeros
362  inline HYPRE_Int NNZ() const { return A->num_nonzeros; }
363  /// Returns the row partitioning
364  /** See @ref hypre_partitioning_descr "here" for a description of the
365  partitioning array. */
366  inline HYPRE_Int *RowPart() { return A->row_starts; }
367  /// Returns the column partitioning
368  /** See @ref hypre_partitioning_descr "here" for a description of the
369  partitioning array. */
370  inline HYPRE_Int *ColPart() { return A->col_starts; }
371  /// Returns the row partitioning (const version)
372  /** See @ref hypre_partitioning_descr "here" for a description of the
373  partitioning array. */
374  inline const HYPRE_Int *RowPart() const { return A->row_starts; }
375  /// Returns the column partitioning (const version)
376  /** See @ref hypre_partitioning_descr "here" for a description of the
377  partitioning array. */
378  inline const HYPRE_Int *ColPart() const { return A->col_starts; }
379  /// Returns the global number of rows
380  inline HYPRE_Int M() const { return A->global_num_rows; }
381  /// Returns the global number of columns
382  inline HYPRE_Int N() const { return A->global_num_cols; }
383 
384  /// Get the local diagonal of the matrix.
385  void GetDiag(Vector &diag) const;
386  /// Get the local diagonal block. NOTE: 'diag' will not own any data.
387  void GetDiag(SparseMatrix &diag) const;
388  /// Get the local off-diagonal block. NOTE: 'offd' will not own any data.
389  void GetOffd(SparseMatrix &offd, HYPRE_Int* &cmap) const;
390 
391  /** Split the matrix into M x N equally sized blocks of parallel matrices.
392  The size of 'blocks' must already be set to M x N. */
393  void GetBlocks(Array2D<HypreParMatrix*> &blocks,
394  bool interleaved_rows = false,
395  bool interleaved_cols = false) const;
396 
397  /// Returns the transpose of *this
398  HypreParMatrix * Transpose() const;
399 
400  /// Returns the number of rows in the diagonal block of the ParCSRMatrix
401  int GetNumRows() const
402  {
403  return internal::to_int(
404  hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(A)));
405  }
406 
407  /// Returns the number of columns in the diagonal block of the ParCSRMatrix
408  int GetNumCols() const
409  {
410  return internal::to_int(
411  hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(A)));
412  }
413 
414  /// Return the global number of rows
415  HYPRE_Int GetGlobalNumRows() const
416  { return hypre_ParCSRMatrixGlobalNumRows(A); }
417 
418  /// Return the global number of columns
419  HYPRE_Int GetGlobalNumCols() const
420  { return hypre_ParCSRMatrixGlobalNumCols(A); }
421 
422  /// Return the parallel row partitioning array.
423  /** See @ref hypre_partitioning_descr "here" for a description of the
424  partitioning array. */
425  HYPRE_Int *GetRowStarts() const { return hypre_ParCSRMatrixRowStarts(A); }
426 
427  /// Return the parallel column partitioning array.
428  /** See @ref hypre_partitioning_descr "here" for a description of the
429  partitioning array. */
430  HYPRE_Int *GetColStarts() const { return hypre_ParCSRMatrixColStarts(A); }
431 
432  /// Computes y = alpha * A * x + beta * y
433  HYPRE_Int Mult(HypreParVector &x, HypreParVector &y,
434  double alpha = 1.0, double beta = 0.0);
435  /// Computes y = alpha * A * x + beta * y
436  HYPRE_Int Mult(HYPRE_ParVector x, HYPRE_ParVector y,
437  double alpha = 1.0, double beta = 0.0);
438  /// Computes y = alpha * A^t * x + beta * y
440  double alpha = 1.0, double beta = 0.0);
441 
442  void Mult(double a, const Vector &x, double b, Vector &y) const;
443  void MultTranspose(double a, const Vector &x, double b, Vector &y) const;
444 
445  virtual void Mult(const Vector &x, Vector &y) const
446  { Mult(1.0, x, 0.0, y); }
447  virtual void MultTranspose(const Vector &x, Vector &y) const
448  { MultTranspose(1.0, x, 0.0, y); }
449 
450  /// Computes y = a * |A| * x + b * y, using entry-wise absolute values of matrix A
451  void AbsMult(double a, const Vector &x, double b, Vector &y) const;
452 
453  /// Computes y = a * |At| * x + b * y, using entry-wise absolute values of the transpose of matrix A
454  void AbsMultTranspose(double a, const Vector &x, double b, Vector &y) const;
455 
456  /** The "Boolean" analog of y = alpha * A * x + beta * y, where elements in
457  the sparsity pattern of the matrix are treated as "true". */
458  void BooleanMult(int alpha, const int *x, int beta, int *y)
459  {
460  internal::hypre_ParCSRMatrixBooleanMatvec(A, alpha, const_cast<int*>(x),
461  beta, y);
462  }
463 
464  /** The "Boolean" analog of y = alpha * A^T * x + beta * y, where elements in
465  the sparsity pattern of the matrix are treated as "true". */
466  void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
467  {
468  internal::hypre_ParCSRMatrixBooleanMatvecT(A, alpha, const_cast<int*>(x),
469  beta, y);
470  }
471 
472  /// Initialize all entries with value.
473  HypreParMatrix &operator=(double value)
474  { internal::hypre_ParCSRMatrixSetConstantValues(A, value); return *this; }
475 
476  /** Perform the operation `*this += B`, assuming that both matrices use the
477  same row and column partitions and the same col_map_offd arrays, or B has
478  an empty off-diagonal block. We also assume that the sparsity pattern of
479  `*this` contains that of `B`. */
480  HypreParMatrix &operator+=(const HypreParMatrix &B) { return Add(1.0, B); }
481 
482  /** Perform the operation `*this += beta*B`, assuming that both matrices use
483  the same row and column partitions and the same col_map_offd arrays, or
484  B has an empty off-diagonal block. We also assume that the sparsity
485  pattern of `*this` contains that of `B`. For a more general case consider
486  the stand-alone function ParAdd described below. */
487  HypreParMatrix &Add(const double beta, const HypreParMatrix &B)
488  {
489  MFEM_VERIFY(internal::hypre_ParCSRMatrixSum(A, beta, B.A) == 0,
490  "error in hypre_ParCSRMatrixSum");
491  return *this;
492  }
493 
494  /** @brief Multiply the HypreParMatrix on the left by a block-diagonal
495  parallel matrix @a D and return the result as a new HypreParMatrix. */
496  /** If @a D has a different number of rows than @a A (this matrix), @a D's
497  row starts array needs to be given (as returned by the methods
498  GetDofOffsets/GetTrueDofOffsets of ParFiniteElementSpace). The new
499  matrix @a D*A uses copies of the row- and column-starts arrays, so "this"
500  matrix and @a row_starts can be deleted.
501  @note This operation is local and does not require communication. */
503  HYPRE_Int* row_starts = NULL) const;
504 
505  /// Scale the local row i by s(i).
506  void ScaleRows(const Vector & s);
507  /// Scale the local row i by 1./s(i)
508  void InvScaleRows(const Vector & s);
509  /// Scale all entries by s: A_scaled = s*A.
510  void operator*=(double s);
511 
512  /// Remove values smaller in absolute value than some threshold
513  void Threshold(double threshold = 0.0);
514 
515  /// If a row contains only zeros, set its diagonal to 1.
516  void EliminateZeroRows() { hypre_ParCSRMatrixFixZeroRows(A); }
517 
518  /** Eliminate rows and columns from the matrix, and rows from the vector B.
519  Modify B with the BC values in X. */
520  void EliminateRowsCols(const Array<int> &rows_cols, const HypreParVector &X,
521  HypreParVector &B);
522 
523  /** Eliminate rows and columns from the matrix and store the eliminated
524  elements in a new matrix Ae (returned), so that the modified matrix and
525  Ae sum to the original matrix. */
526  HypreParMatrix* EliminateRowsCols(const Array<int> &rows_cols);
527 
528  /** Eliminate columns from the matrix and store the eliminated elements in a
529  new matrix Ae (returned) so that the modified matrix and Ae sum to the
530  original matrix. */
532 
533  /// Eliminate rows from the diagonal and off-diagonal blocks of the matrix.
534  void EliminateRows(const Array<int> &rows);
535 
536  /// Prints the locally owned rows in parallel
537  void Print(const char *fname, HYPRE_Int offi = 0, HYPRE_Int offj = 0);
538  /// Reads the matrix from a file
539  void Read(MPI_Comm comm, const char *fname);
540  /// Read a matrix saved as a HYPRE_IJMatrix
541  void Read_IJMatrix(MPI_Comm comm, const char *fname);
542 
543  /// Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
544  void PrintCommPkg(std::ostream &out = mfem::out) const;
545 
546  /// Calls hypre's destroy function
547  virtual ~HypreParMatrix() { Destroy(); }
548 
549  Type GetType() const { return Hypre_ParCSR; }
550 };
551 
552 /** @brief Return a new matrix `C = alpha*A + beta*B`, assuming that both `A`
553  and `B` use the same row and column partitions and the same `col_map_offd`
554  arrays. */
555 HypreParMatrix *Add(double alpha, const HypreParMatrix &A,
556  double beta, const HypreParMatrix &B);
557 
558 /** Returns the matrix @a A * @a B. Returned matrix does not necessarily own
559  row or column starts unless the bool @a own_matrix is set to true. */
560 HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B,
561  bool own_matrix = false);
562 /// Returns the matrix A + B
563 /** It is assumed that both matrices use the same row and column partitions and
564  the same col_map_offd arrays. */
565 HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B);
566 
567 /// Returns the matrix P^t * A * P
568 HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P);
569 /// Returns the matrix Rt^t * A * P
570 HypreParMatrix * RAP(const HypreParMatrix * Rt, const HypreParMatrix *A,
571  const HypreParMatrix *P);
572 
573 /// Returns a merged hypre matrix constructed from hypre matrix blocks.
574 /** It is assumed that all block matrices use the same communicator, and the
575  block sizes are consistent in rows and columns. Rows and columns are
576  renumbered but not redistributed in parallel, e.g. the block rows owned by
577  each process remain on that process in the resulting matrix. Some blocks can
578  be NULL. Each block and the entire system can be rectangular. Scalability to
579  extremely large processor counts is limited by global MPI communication, see
580  GatherBlockOffsetData() in hypre.cpp. */
581 HypreParMatrix * HypreParMatrixFromBlocks(Array2D<HypreParMatrix*> &blocks,
582  Array2D<double> *blockCoeff=NULL);
583 
584 /** Eliminate essential BC specified by 'ess_dof_list' from the solution X to
585  the r.h.s. B. Here A is a matrix with eliminated BC, while Ae is such that
586  (A+Ae) is the original (Neumann) matrix before elimination. */
587 void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae,
588  const Array<int> &ess_dof_list, const Vector &X, Vector &B);
589 
590 
591 /// Parallel smoothers in hypre
592 class HypreSmoother : public Solver
593 {
594 protected:
595  /// The linear system matrix
597  /// Right-hand side and solution vectors
598  mutable HypreParVector *B, *X;
599  /// Temporary vectors
600  mutable HypreParVector *V, *Z;
601  /// FIR Filter Temporary Vectors
602  mutable HypreParVector *X0, *X1;
603 
604  /** Smoother type from hypre_ParCSRRelax() in ams.c plus extensions, see the
605  enumeration Type below. */
606  int type;
607  /// Number of relaxation sweeps
609  /// Damping coefficient (usually <= 1)
610  double relax_weight;
611  /// SOR parameter (usually in (0,2))
612  double omega;
613  /// Order of the smoothing polynomial
615  /// Fraction of spectrum to smooth for polynomial relaxation
617  /// Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}
619 
620  /// Taubin's lambda-mu method parameters
621  double lambda;
622  double mu;
624 
625  /// l1 norms of the rows of A
626  double *l1_norms;
627  /// If set, take absolute values of the computed l1_norms
629  /// Number of CG iterations to determine eigenvalue estimates
631  /// Maximal eigenvalue estimate for polynomial smoothing
632  double max_eig_est;
633  /// Minimal eigenvalue estimate for polynomial smoothing
634  double min_eig_est;
635  /// Parameters for windowing function of FIR filter
636  double window_params[3];
637 
638  /// Combined coefficients for windowing and Chebyshev polynomials.
639  double* fir_coeffs;
640 
641 public:
642  /** Hypre smoother types:
643  0 = Jacobi
644  1 = l1-scaled Jacobi
645  2 = l1-scaled block Gauss-Seidel/SSOR
646  4 = truncated l1-scaled block Gauss-Seidel/SSOR
647  5 = lumped Jacobi
648  6 = Gauss-Seidel
649  16 = Chebyshev
650  1001 = Taubin polynomial smoother
651  1002 = FIR polynomial smoother. */
652  enum Type { Jacobi = 0, l1Jacobi = 1, l1GS = 2, l1GStr = 4, lumpedJacobi = 5,
653  GS = 6, Chebyshev = 16, Taubin = 1001, FIR = 1002
654  };
655 
656  HypreSmoother();
657 
659  int relax_times = 1, double relax_weight = 1.0,
660  double omega = 1.0, int poly_order = 2,
661  double poly_fraction = .3, int eig_est_cg_iter = 10);
662 
663  /// Set the relaxation type and number of sweeps
665  /// Set SOR-related parameters
666  void SetSOROptions(double relax_weight, double omega);
667  /// Set parameters for polynomial smoothing
668  /** By default, 10 iterations of CG are used to estimate the eigenvalues.
669  Setting eig_est_cg_iter = 0 uses hypre's hypre_ParCSRMaxEigEstimate() instead. */
670  void SetPolyOptions(int poly_order, double poly_fraction,
671  int eig_est_cg_iter = 10);
672  /// Set parameters for Taubin's lambda-mu method
673  void SetTaubinOptions(double lambda, double mu, int iter);
674 
675  /// Convenience function for setting canonical windowing parameters
676  void SetWindowByName(const char* window_name);
677  /// Set parameters for windowing function for FIR smoother.
678  void SetWindowParameters(double a, double b, double c);
679  /// Compute window and Chebyshev coefficients for given polynomial order.
680  void SetFIRCoefficients(double max_eig);
681 
682  /// After computing l1-norms, replace them with their absolute values.
683  /** By default, the l1-norms take their sign from the corresponding diagonal
684  entries in the associated matrix. */
685  void SetPositiveDiagonal(bool pos = true) { pos_l1_norms = pos; }
686 
687  /** Set/update the associated operator. Must be called after setting the
688  HypreSmoother type and options. */
689  virtual void SetOperator(const Operator &op);
690 
691  /// Relax the linear system Ax=b
692  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
693  virtual void Mult(const Vector &b, Vector &x) const;
694 
695  virtual ~HypreSmoother();
696 };
697 
698 
699 /// Abstract class for hypre's solvers and preconditioners
700 class HypreSolver : public Solver
701 {
702 public:
703  /// How to treat errors returned by hypre function calls.
705  {
706  IGNORE_HYPRE_ERRORS, ///< Ignore hypre errors (see e.g. HypreADS)
707  WARN_HYPRE_ERRORS, ///< Issue warnings on hypre errors
708  ABORT_HYPRE_ERRORS ///< Abort on hypre errors (default in base class)
709  };
710 
711 protected:
712  /// The linear system matrix
714 
715  /// Right-hand side and solution vector
716  mutable HypreParVector *B, *X;
717 
718  /// Was hypre's Setup function called already?
719  mutable int setup_called;
720 
721  /// How to treat hypre errors.
723 
724 public:
725  HypreSolver();
726 
728 
729  /// Typecast to HYPRE_Solver -- return the solver
730  virtual operator HYPRE_Solver() const = 0;
731 
732  /// hypre's internal Setup function
733  virtual HYPRE_PtrToParSolverFcn SetupFcn() const = 0;
734  /// hypre's internal Solve function
735  virtual HYPRE_PtrToParSolverFcn SolveFcn() const = 0;
736 
737  virtual void SetOperator(const Operator &op)
738  { mfem_error("HypreSolvers do not support SetOperator!"); }
739 
740  /// Solve the linear system Ax=b
741  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
742  virtual void Mult(const Vector &b, Vector &x) const;
743 
744  /** @brief Set the behavior for treating hypre errors, see the ErrorMode
745  enum. The default mode in the base class is ABORT_HYPRE_ERRORS. */
746  /** Currently, there are three cases in derived classes where the error flag
747  is set to IGNORE_HYPRE_ERRORS:
748  * in the method HypreBoomerAMG::SetElasticityOptions(), and
749  * in the constructor of classes HypreAMS and HypreADS.
750  The reason for this is that a nonzero hypre error is returned) when
751  hypre_ParCSRComputeL1Norms() encounters zero row in a matrix, which is
752  expected in some cases with the above solvers. */
753  void SetErrorMode(ErrorMode err_mode) const { error_mode = err_mode; }
754 
755  virtual ~HypreSolver();
756 };
757 
758 /// PCG solver in hypre
759 class HyprePCG : public HypreSolver
760 {
761 private:
762  HYPRE_Solver pcg_solver;
763 
764  HypreSolver * precond;
765 
766 public:
767  HyprePCG(MPI_Comm comm);
768 
770 
771  virtual void SetOperator(const Operator &op);
772 
773  void SetTol(double tol);
774  void SetMaxIter(int max_iter);
775  void SetLogging(int logging);
776  void SetPrintLevel(int print_lvl);
777 
778  /// Set the hypre solver to be used as a preconditioner
779  void SetPreconditioner(HypreSolver &precond);
780 
781  /** Use the L2 norm of the residual for measuring PCG convergence, plus
782  (optionally) 1) periodically recompute true residuals from scratch; and
783  2) enable residual-based stopping criteria. */
784  void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0);
785 
786  /// deprecated: use SetZeroInitialIterate()
787  MFEM_DEPRECATED void SetZeroInintialIterate() { iterative_mode = false; }
788 
789  /// non-hypre setting
791 
792  void GetNumIterations(int &num_iterations)
793  {
794  HYPRE_Int num_it;
795  HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_it);
796  num_iterations = internal::to_int(num_it);
797  }
798 
799  /// The typecast to HYPRE_Solver returns the internal pcg_solver
800  virtual operator HYPRE_Solver() const { return pcg_solver; }
801 
802  /// PCG Setup function
803  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
804  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSetup; }
805  /// PCG Solve function
806  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
807  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSolve; }
808 
809  /// Solve Ax=b with hypre's PCG
810  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
811  using HypreSolver::Mult;
812 
813  virtual ~HyprePCG();
814 };
815 
816 /// GMRES solver in hypre
817 class HypreGMRES : public HypreSolver
818 {
819 private:
820  HYPRE_Solver gmres_solver;
821 
822  HypreSolver * precond;
823 
824  /// Default, generally robust, GMRES options
825  void SetDefaultOptions();
826 
827 public:
828  HypreGMRES(MPI_Comm comm);
829 
831 
832  virtual void SetOperator(const Operator &op);
833 
834  void SetTol(double tol);
835  void SetMaxIter(int max_iter);
836  void SetKDim(int dim);
837  void SetLogging(int logging);
838  void SetPrintLevel(int print_lvl);
839 
840  /// Set the hypre solver to be used as a preconditioner
841  void SetPreconditioner(HypreSolver &precond);
842 
843  /// deprecated: use SetZeroInitialIterate()
844  MFEM_DEPRECATED void SetZeroInintialIterate() { iterative_mode = false; }
845 
846  /// non-hypre setting
848 
849  /// The typecast to HYPRE_Solver returns the internal gmres_solver
850  virtual operator HYPRE_Solver() const { return gmres_solver; }
851 
852  /// GMRES Setup function
853  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
854  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSetup; }
855  /// GMRES Solve function
856  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
857  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSolve; }
858 
859  /// Solve Ax=b with hypre's GMRES
860  virtual void Mult (const HypreParVector &b, HypreParVector &x) const;
861  using HypreSolver::Mult;
862 
863  virtual ~HypreGMRES();
864 };
865 
866 /// Flexible GMRES solver in hypre
867 class HypreFGMRES : public HypreSolver
868 {
869 private:
870  HYPRE_Solver fgmres_solver;
871 
872  HypreSolver * precond;
873 
874  /// Default, generally robust, FGMRES options
875  void SetDefaultOptions();
876 
877 public:
878  HypreFGMRES(MPI_Comm comm);
879 
881 
882  virtual void SetOperator(const Operator &op);
883 
884  void SetTol(double tol);
885  void SetMaxIter(int max_iter);
886  void SetKDim(int dim);
887  void SetLogging(int logging);
888  void SetPrintLevel(int print_lvl);
889 
890  /// Set the hypre solver to be used as a preconditioner
891  void SetPreconditioner(HypreSolver &precond);
892 
893  /// deprecated: use SetZeroInitialIterate()
894  MFEM_DEPRECATED void SetZeroInintialIterate() { iterative_mode = false; }
895 
896  /// non-hypre setting
898 
899  /// The typecast to HYPRE_Solver returns the internal fgmres_solver
900  virtual operator HYPRE_Solver() const { return fgmres_solver; }
901 
902  /// FGMRES Setup function
903  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
904  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRFlexGMRESSetup; }
905  /// FGMRES Solve function
906  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
907  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRFlexGMRESSolve; }
908 
909  /// Solve Ax=b with hypre's FGMRES
910  virtual void Mult (const HypreParVector &b, HypreParVector &x) const;
911  using HypreSolver::Mult;
912 
913  virtual ~HypreFGMRES();
914 };
915 
916 /// The identity operator as a hypre solver
918 {
919 public:
920  virtual operator HYPRE_Solver() const { return NULL; }
921 
922  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
923  { return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentitySetup; }
924  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
925  { return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentity; }
926 
927  virtual ~HypreIdentity() { }
928 };
929 
930 /// Jacobi preconditioner in hypre
932 {
933 public:
936  virtual operator HYPRE_Solver() const { return NULL; }
937 
938  virtual void SetOperator(const Operator &op);
939 
940  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
941  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScaleSetup; }
942  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
943  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScale; }
944 
945  HypreParMatrix* GetData() { return A; }
946  virtual ~HypreDiagScale() { }
947 };
948 
949 /// The ParaSails preconditioner in hypre
951 {
952 private:
953  HYPRE_Solver sai_precond;
954 
955  /// Default, generally robust, ParaSails options
956  void SetDefaultOptions();
957 
958  // If sai_precond is NULL, this method allocates it and sets default options.
959  // Otherwise the method saves the options from sai_precond, destroys it,
960  // allocates a new object, and sets its options to the saved values.
961  void ResetSAIPrecond(MPI_Comm comm);
962 
963 public:
964  HypreParaSails(MPI_Comm comm);
965 
967 
968  virtual void SetOperator(const Operator &op);
969 
970  void SetSymmetry(int sym);
971 
972  /// The typecast to HYPRE_Solver returns the internal sai_precond
973  virtual operator HYPRE_Solver() const { return sai_precond; }
974 
975  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
976  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSetup; }
977  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
978  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSolve; }
979 
980  virtual ~HypreParaSails();
981 };
982 
983 /** The Euclid preconditioner in Hypre
984 
985  Euclid implements the Parallel Incomplete LU factorization technique. For
986  more information see:
987 
988  "A Scalable Parallel Algorithm for Incomplete Factor Preconditioning" by
989  David Hysom and Alex Pothen, https://doi.org/10.1137/S1064827500376193
990 */
991 class HypreEuclid : public HypreSolver
992 {
993 private:
994  HYPRE_Solver euc_precond;
995 
996  /// Default, generally robust, Euclid options
997  void SetDefaultOptions();
998 
999  // If euc_precond is NULL, this method allocates it and sets default options.
1000  // Otherwise the method saves the options from euc_precond, destroys it,
1001  // allocates a new object, and sets its options to the saved values.
1002  void ResetEuclidPrecond(MPI_Comm comm);
1003 
1004 public:
1005  HypreEuclid(MPI_Comm comm);
1006 
1008 
1009  virtual void SetOperator(const Operator &op);
1010 
1011  /// The typecast to HYPRE_Solver returns the internal euc_precond
1012  virtual operator HYPRE_Solver() const { return euc_precond; }
1013 
1014  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1015  { return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSetup; }
1016  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1017  { return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSolve; }
1018 
1019  virtual ~HypreEuclid();
1020 };
1021 
1022 #if MFEM_HYPRE_VERSION >= 21900
1023 /**
1024 @brief Wrapper for Hypre's native parallel ILU preconditioner.
1025 
1026 The default ILU factorization type is ILU(k). If you need to change this, or
1027 any other option, you can use the HYPRE_Solver method to cast the object for use
1028 with Hypre's native functions. For example, if want to use natural ordering
1029 rather than RCM reordering, you can use the following approach:
1030 
1031 @code
1032 mfem::HypreILU ilu();
1033 int reorder_type = 0;
1034 HYPRE_ILUSetLocalReordering(ilu, reorder_type);
1035 @endcode
1036 */
1037 class HypreILU : public HypreSolver
1038 {
1039 private:
1040  HYPRE_Solver ilu_precond;
1041 
1042  /// Set the ILU default options
1043  void SetDefaultOptions();
1044 
1045  /** Reset the ILU preconditioner.
1046  @note If ilu_precond is NULL, this method allocates; otherwise it destroys
1047  ilu_precond and allocates a new object. In both cases the default options
1048  are set. */
1049  void ResetILUPrecond();
1050 
1051 public:
1052  /// Constructor; sets the default options
1053  HypreILU();
1054 
1055  virtual ~HypreILU();
1056 
1057  /// Set the fill level for ILU(k); the default is k=1.
1058  void SetLevelOfFill(HYPRE_Int lev_fill);
1059 
1060  /// Set the print level: 0 = none, 1 = setup, 2 = solve, 3 = setup+solve
1061  void SetPrintLevel(HYPRE_Int print_level);
1062 
1063  /// The typecast to HYPRE_Solver returns the internal ilu_precond
1064  virtual operator HYPRE_Solver() const { return ilu_precond; }
1065 
1066  virtual void SetOperator(const Operator &op);
1067 
1068  /// ILU Setup function
1069  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1070  { return (HYPRE_PtrToParSolverFcn) HYPRE_ILUSetup; }
1071 
1072  /// ILU Solve function
1073  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1074  { return (HYPRE_PtrToParSolverFcn) HYPRE_ILUSolve; }
1075 };
1076 #endif
1077 
1078 /// The BoomerAMG solver in hypre
1080 {
1081 private:
1082  HYPRE_Solver amg_precond;
1083 
1084  /// Rigid body modes
1086 
1087  /// Finite element space for elasticity problems, see SetElasticityOptions()
1088  ParFiniteElementSpace *fespace;
1089 
1090  /// Recompute the rigid-body modes vectors (in the rbms array)
1091  void RecomputeRBMs();
1092 
1093  /// Default, generally robust, BoomerAMG options
1094  void SetDefaultOptions();
1095 
1096  // If amg_precond is NULL, allocates it and sets default options.
1097  // Otherwise saves the options from amg_precond, destroys it, allocates a new
1098  // one, and sets its options to the saved values.
1099  void ResetAMGPrecond();
1100 
1101 public:
1102  HypreBoomerAMG();
1103 
1105 
1106  virtual void SetOperator(const Operator &op);
1107 
1108  /** More robust options for systems, such as elasticity. */
1109  void SetSystemsOptions(int dim, bool order_bynodes=false);
1110 
1111  /** A special elasticity version of BoomerAMG that takes advantage of
1112  geometric rigid body modes and could perform better on some problems, see
1113  "Improving algebraic multigrid interpolation operators for linear
1114  elasticity problems", Baker, Kolev, Yang, NLAA 2009, DOI:10.1002/nla.688.
1115  This solver assumes Ordering::byVDIM in the FiniteElementSpace used to
1116  construct A. */
1118 
1119  void SetPrintLevel(int print_level)
1120  { HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level); }
1121 
1122  /// The typecast to HYPRE_Solver returns the internal amg_precond
1123  virtual operator HYPRE_Solver() const { return amg_precond; }
1124 
1125  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1126  { return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSetup; }
1127  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1128  { return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSolve; }
1129 
1130  virtual ~HypreBoomerAMG();
1131 };
1132 
1133 /// Compute the discrete gradient matrix between the nodal linear and ND1 spaces
1134 HypreParMatrix* DiscreteGrad(ParFiniteElementSpace *edge_fespace,
1135  ParFiniteElementSpace *vert_fespace);
1136 /// Compute the discrete curl matrix between the ND1 and RT0 spaces
1137 HypreParMatrix* DiscreteCurl(ParFiniteElementSpace *face_fespace,
1138  ParFiniteElementSpace *edge_fespace);
1139 
1140 /// The Auxiliary-space Maxwell Solver in hypre
1141 class HypreAMS : public HypreSolver
1142 {
1143 private:
1144  /// Constuct AMS solver from finite element space
1145  void Init(ParFiniteElementSpace *edge_space);
1146 
1147  HYPRE_Solver ams;
1148 
1149  /// Vertex coordinates
1150  HypreParVector *x, *y, *z;
1151  /// Discrete gradient matrix
1152  HypreParMatrix *G;
1153  /// Nedelec interpolation matrix and its components
1154  HypreParMatrix *Pi, *Pix, *Piy, *Piz;
1155 
1156 public:
1157  HypreAMS(ParFiniteElementSpace *edge_fespace);
1158 
1160 
1161  virtual void SetOperator(const Operator &op);
1162 
1163  void SetPrintLevel(int print_lvl);
1164 
1165  /// Set this option when solving a curl-curl problem with zero mass term
1166  void SetSingularProblem() { HYPRE_AMSSetBetaPoissonMatrix(ams, NULL); }
1167 
1168  /// The typecast to HYPRE_Solver returns the internal ams object
1169  virtual operator HYPRE_Solver() const { return ams; }
1170 
1171  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1172  { return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSetup; }
1173  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1174  { return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSolve; }
1175 
1176  virtual ~HypreAMS();
1177 };
1178 
1179 /// The Auxiliary-space Divergence Solver in hypre
1180 class HypreADS : public HypreSolver
1181 {
1182 private:
1183  /// Constuct ADS solver from finite element space
1184  void Init(ParFiniteElementSpace *face_fespace);
1185 
1186  HYPRE_Solver ads;
1187 
1188  /// Vertex coordinates
1189  HypreParVector *x, *y, *z;
1190  /// Discrete gradient matrix
1191  HypreParMatrix *G;
1192  /// Discrete curl matrix
1193  HypreParMatrix *C;
1194  /// Nedelec interpolation matrix and its components
1195  HypreParMatrix *ND_Pi, *ND_Pix, *ND_Piy, *ND_Piz;
1196  /// Raviart-Thomas interpolation matrix and its components
1197  HypreParMatrix *RT_Pi, *RT_Pix, *RT_Piy, *RT_Piz;
1198 
1199 public:
1200  HypreADS(ParFiniteElementSpace *face_fespace);
1201 
1203 
1204  virtual void SetOperator(const Operator &op);
1205 
1206  void SetPrintLevel(int print_lvl);
1207 
1208  /// The typecast to HYPRE_Solver returns the internal ads object
1209  virtual operator HYPRE_Solver() const { return ads; }
1210 
1211  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1212  { return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSetup; }
1213  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1214  { return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSolve; }
1215 
1216  virtual ~HypreADS();
1217 };
1218 
1219 /** LOBPCG eigenvalue solver in hypre
1220 
1221  The Locally Optimal Block Preconditioned Conjugate Gradient (LOBPCG)
1222  eigenvalue solver is designed to find the lowest eigenmodes of the
1223  generalized eigenvalue problem:
1224  A x = lambda M x
1225  where A is symmetric, potentially indefinite and M is symmetric positive
1226  definite. The eigenvectors are M-orthonormal, meaning that
1227  x^T M x = 1 and x^T M y = 0,
1228  if x and y are distinct eigenvectors. The matrix M is optional and is
1229  assumed to be the identity if left unset.
1230 
1231  The efficiency of LOBPCG relies on the availability of a suitable
1232  preconditioner for the matrix A. The preconditioner is supplied through the
1233  SetPreconditioner() method. It should be noted that the operator used with
1234  the preconditioner need not be A itself.
1235 
1236  For more information regarding LOBPCG see "Block Locally Optimal
1237  Preconditioned Eigenvalue Xolvers (BLOPEX) in Hypre and PETSc" by
1238  A. Knyazev, M. Argentati, I. Lashuk, and E. Ovtchinnikov, SISC, 29(5),
1239  2224-2239, 2007.
1240 */
1242 {
1243 private:
1244  MPI_Comm comm;
1245  int myid;
1246  int numProcs;
1247  int nev; // Number of desired eigenmodes
1248  int seed; // Random seed used for initial vectors
1249 
1250  HYPRE_Int glbSize; // Global number of DoFs in the linear system
1251  HYPRE_Int * part; // Row partitioning of the linear system
1252 
1253  // Pointer to HYPRE's solver struct
1254  HYPRE_Solver lobpcg_solver;
1255 
1256  // Interface for matrix storage type
1257  mv_InterfaceInterpreter interpreter;
1258 
1259  // Interface for setting up and performing matrix-vector products
1260  HYPRE_MatvecFunctions matvec_fn;
1261 
1262  // Eigenvalues
1263  Array<double> eigenvalues;
1264 
1265  // Forward declaration
1266  class HypreMultiVector;
1267 
1268  // MultiVector to store eigenvectors
1269  HypreMultiVector * multi_vec;
1270 
1271  // Empty vectors used to setup the matrices and preconditioner
1272  HypreParVector * x;
1273 
1274  // An optional operator which projects vectors into a desired subspace
1275  Operator * subSpaceProj;
1276 
1277  /// Internal class to represent a set of eigenvectors
1278  class HypreMultiVector
1279  {
1280  private:
1281  // Pointer to hypre's multi-vector object
1282  mv_MultiVectorPtr mv_ptr;
1283 
1284  // Wrappers for each member of the multivector
1285  HypreParVector ** hpv;
1286 
1287  // Number of vectors in the multivector
1288  int nv;
1289 
1290  public:
1291  HypreMultiVector(int n, HypreParVector & v,
1292  mv_InterfaceInterpreter & interpreter);
1293  ~HypreMultiVector();
1294 
1295  /// Set random values
1296  void Randomize(HYPRE_Int seed);
1297 
1298  /// Extract a single HypreParVector object
1299  HypreParVector & GetVector(unsigned int i);
1300 
1301  /// Transfers ownership of data to returned array of vectors
1302  HypreParVector ** StealVectors();
1303 
1304  operator mv_MultiVectorPtr() const { return mv_ptr; }
1305 
1306  mv_MultiVectorPtr & GetMultiVector() { return mv_ptr; }
1307  };
1308 
1309  static void * OperatorMatvecCreate( void *A, void *x );
1310  static HYPRE_Int OperatorMatvec( void *matvec_data,
1311  HYPRE_Complex alpha,
1312  void *A,
1313  void *x,
1314  HYPRE_Complex beta,
1315  void *y );
1316  static HYPRE_Int OperatorMatvecDestroy( void *matvec_data );
1317 
1318  static HYPRE_Int PrecondSolve(void *solver,
1319  void *A,
1320  void *b,
1321  void *x);
1322  static HYPRE_Int PrecondSetup(void *solver,
1323  void *A,
1324  void *b,
1325  void *x);
1326 
1327 public:
1328  HypreLOBPCG(MPI_Comm comm);
1329  ~HypreLOBPCG();
1330 
1331  void SetTol(double tol);
1332  void SetRelTol(double rel_tol);
1333  void SetMaxIter(int max_iter);
1334  void SetPrintLevel(int logging);
1335  void SetNumModes(int num_eigs) { nev = num_eigs; }
1336  void SetPrecondUsageMode(int pcg_mode);
1337  void SetRandomSeed(int s) { seed = s; }
1338  void SetInitialVectors(int num_vecs, HypreParVector ** vecs);
1339 
1340  // The following four methods support general operators
1341  void SetPreconditioner(Solver & precond);
1342  void SetOperator(Operator & A);
1343  void SetMassMatrix(Operator & M);
1344  void SetSubSpaceProjector(Operator & proj) { subSpaceProj = &proj; }
1345 
1346  /// Solve the eigenproblem
1347  void Solve();
1348 
1349  /// Collect the converged eigenvalues
1350  void GetEigenvalues(Array<double> & eigenvalues);
1351 
1352  /// Extract a single eigenvector
1353  HypreParVector & GetEigenvector(unsigned int i);
1354 
1355  /// Transfer ownership of the converged eigenvectors
1356  HypreParVector ** StealEigenvectors() { return multi_vec->StealVectors(); }
1357 };
1358 
1359 /** AME eigenvalue solver in hypre
1360 
1361  The Auxiliary space Maxwell Eigensolver (AME) is designed to find
1362  the lowest eigenmodes of the generalized eigenvalue problem:
1363  Curl Curl x = lambda M x
1364  where the Curl Curl operator is discretized using Nedelec finite element
1365  basis functions. Properties of this discretization are essential to
1366  eliminating the large null space of the Curl Curl operator.
1367 
1368  This eigensolver relies upon the LOBPCG eigensolver internally. It is also
1369  expected that the preconditioner supplied to this method will be the
1370  HypreAMS preconditioner defined above.
1371 
1372  As with LOBPCG, the operator set in the preconditioner need not be the same
1373  as A. This flexibility may be useful in solving eigenproblems which bare a
1374  strong resemblance to the Curl Curl problems for which AME is designed.
1375 
1376  Unlike LOBPCG, this eigensolver requires that the mass matrix be set.
1377  It is possible to circumvent this by passing an identity operator as the
1378  mass matrix but it seems unlikely that this would be useful so it is not the
1379  default behavior.
1380 */
1382 {
1383 private:
1384  int myid;
1385  int numProcs;
1386  int nev; // Number of desired eigenmodes
1387  bool setT;
1388 
1389  // Pointer to HYPRE's AME solver struct
1390  HYPRE_Solver ame_solver;
1391 
1392  // Pointer to HYPRE's AMS solver struct
1393  HypreSolver * ams_precond;
1394 
1395  // Eigenvalues
1396  HYPRE_Real * eigenvalues;
1397 
1398  // MultiVector to store eigenvectors
1399  HYPRE_ParVector * multi_vec;
1400 
1401  HypreParVector ** eigenvectors;
1402 
1403  void createDummyVectors();
1404 
1405 public:
1406  HypreAME(MPI_Comm comm);
1407  ~HypreAME();
1408 
1409  void SetTol(double tol);
1410  void SetRelTol(double rel_tol);
1411  void SetMaxIter(int max_iter);
1412  void SetPrintLevel(int logging);
1413  void SetNumModes(int num_eigs);
1414 
1415  // The following four methods support operators of type HypreParMatrix.
1416  void SetPreconditioner(HypreSolver & precond);
1417  void SetOperator(HypreParMatrix & A);
1418  void SetMassMatrix(HypreParMatrix & M);
1419 
1420  /// Solve the eigenproblem
1421  void Solve();
1422 
1423  /// Collect the converged eigenvalues
1424  void GetEigenvalues(Array<double> & eigenvalues);
1425 
1426  /// Extract a single eigenvector
1427  HypreParVector & GetEigenvector(unsigned int i);
1428 
1429  /// Transfer ownership of the converged eigenvectors
1431 };
1432 
1433 }
1434 
1435 #endif // MFEM_USE_MPI
1436 
1437 #endif
Type GetType() const
Definition: hypre.hpp:549
virtual ~HypreBoomerAMG()
Definition: hypre.cpp:3603
const HYPRE_Int * ColPart() const
Returns the column partitioning (const version)
Definition: hypre.hpp:378
void SetTol(double tol)
Definition: hypre.cpp:2619
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2820
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:1368
HypreADS(ParFiniteElementSpace *face_fespace)
Definition: hypre.cpp:3844
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:330
HypreEuclid(MPI_Comm comm)
Definition: hypre.cpp:3176
void SetErrorMode(ErrorMode err_mode) const
Set the behavior for treating hypre errors, see the ErrorMode enum. The default mode in the base clas...
Definition: hypre.hpp:753
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2935
Vector * GlobalVector() const
Returns the global vector in each processor.
Definition: hypre.cpp:131
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1141
HypreParVector * X0
FIR Filter Temporary Vectors.
Definition: hypre.hpp:602
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:4335
double min_eig_est
Minimal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:634
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:924
ErrorMode
How to treat errors returned by hypre function calls.
Definition: hypre.hpp:704
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:1180
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
FGMRES Solve function.
Definition: hypre.hpp:906
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:4116
void Print(const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0)
Prints the locally owned rows in parallel.
Definition: hypre.cpp:1414
MPI_Comm GetComm()
MPI communicator.
Definition: hypre.hpp:114
HYPRE_Int N() const
Returns the global number of columns.
Definition: hypre.hpp:382
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1640
virtual ~HypreEuclid()
Definition: hypre.cpp:3239
int setup_called
Was hypre&#39;s Setup function called already?
Definition: hypre.hpp:719
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2599
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2980
HYPRE_Int * ColPart()
Returns the column partitioning.
Definition: hypre.hpp:370
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to &#39;master&#39;.
Definition: hypre.cpp:785
void Read_IJMatrix(MPI_Comm comm, const char *fname)
Read a matrix saved as a HYPRE_IJMatrix.
Definition: hypre.cpp:1434
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_Int *row_starts=NULL) const
Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result a...
Definition: hypre.cpp:1081
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A^t * x + beta * y.
Definition: hypre.cpp:1045
HypreParMatrix & Add(const double beta, const HypreParMatrix &B)
Definition: hypre.hpp:487
HypreParVector * B
Right-hand side and solution vector.
Definition: hypre.hpp:716
Wrapper for Hypre&#39;s native parallel ILU preconditioner.
Definition: hypre.hpp:1037
HypreDiagScale(HypreParMatrix &A)
Definition: hypre.hpp:935
double window_params[3]
Parameters for windowing function of FIR filter.
Definition: hypre.hpp:636
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: hypre.hpp:445
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1173
void SetInitialVectors(int num_vecs, HypreParVector **vecs)
Definition: hypre.cpp:4341
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2775
void SetWindowByName(const char *window_name)
Convenience function for setting canonical windowing parameters.
Definition: hypre.cpp:2215
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:4602
Issue warnings on hypre errors.
Definition: hypre.hpp:707
void SetPreconditioner(HypreSolver &precond)
Definition: hypre.cpp:4568
HyprePCG(MPI_Comm comm)
Definition: hypre.cpp:2581
void SetTol(double tol)
Definition: hypre.cpp:4537
HypreILU()
Constructor; sets the default options.
Definition: hypre.cpp:3246
HypreAMS(ParFiniteElementSpace *edge_fespace)
Definition: hypre.cpp:3613
void SetKDim(int dim)
Definition: hypre.cpp:2965
virtual ~HypreGMRES()
Definition: hypre.cpp:2898
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2975
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:4078
void SetPreconditioner(Solver &precond)
Definition: hypre.cpp:4258
MFEM_DEPRECATED void SetZeroInintialIterate()
deprecated: use SetZeroInitialIterate()
Definition: hypre.hpp:844
void SetMassMatrix(Operator &M)
Definition: hypre.cpp:4313
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetPrintLevel(int logging)
Definition: hypre.cpp:4243
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:638
virtual ~HypreAMS()
Definition: hypre.cpp:3824
void SetTol(double tol)
Definition: hypre.cpp:2795
virtual ~HypreILU()
Definition: hypre.cpp:3316
void SetOperator(HypreParMatrix &A)
Definition: hypre.cpp:4574
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
GMRES Solve function.
Definition: hypre.hpp:856
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.hpp:1356
HypreParMatrix & operator+=(const HypreParMatrix &B)
Definition: hypre.hpp:480
bool pos_l1_norms
If set, take absolute values of the computed l1_norms.
Definition: hypre.hpp:628
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
Definition: hypre.cpp:1623
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
Definition: hypre.cpp:969
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2815
void SetZeroInitialIterate()
non-hypre setting
Definition: hypre.hpp:897
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
Definition: hypre.cpp:1178
virtual MFEM_DEPRECATED N_Vector ToNVector()
(DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_PARALLEL.
Definition: hypre.cpp:186
double poly_fraction
Fraction of spectrum to smooth for polynomial relaxation.
Definition: hypre.hpp:616
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2634
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:3299
HypreParMatrix * DiscreteGrad(ParFiniteElementSpace *edge_fespace, ParFiniteElementSpace *vert_fespace)
Compute the discrete gradient matrix between the nodal linear and ND1 spaces.
void CopyRowStarts()
Definition: hypre.cpp:809
int poly_scale
Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}.
Definition: hypre.hpp:618
void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
Definition: hypre.hpp:466
void SetTol(double tol)
Definition: hypre.cpp:4221
HypreParVector * X1
Definition: hypre.hpp:602
HypreFGMRES(MPI_Comm comm)
Definition: hypre.cpp:2904
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s FGMRES.
Definition: hypre.cpp:2989
The identity operator as a hypre solver.
Definition: hypre.hpp:917
void SetWindowParameters(double a, double b, double c)
Set parameters for windowing function for FIR smoother.
Definition: hypre.cpp:2230
void EliminateZeroRows()
If a row contains only zeros, set its diagonal to 1.
Definition: hypre.hpp:516
HypreGMRES(MPI_Comm comm)
Definition: hypre.cpp:2744
HYPRE_Int NNZ() const
Returns the global number of nonzeros.
Definition: hypre.hpp:362
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:3143
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s GMRES.
Definition: hypre.cpp:2830
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
Definition: hypre.hpp:685
void Print(const char *fname) const
Prints the locally owned rows in parallel.
Definition: hypre.cpp:171
void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
Definition: hypre.cpp:1981
void AbsMultTranspose(double a, const Vector &x, double b, Vector &y) const
Computes y = a * |At| * x + b * y, using entry-wise absolute values of the transpose of matrix A...
Definition: hypre.cpp:1066
HypreLOBPCG(MPI_Comm comm)
Definition: hypre.cpp:4191
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:1928
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
Definition: hypre.cpp:2369
hypre_ParVector * StealParVector()
Changes the ownership of the the vector.
Definition: hypre.hpp:131
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1079
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:3063
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:415
void SetSymmetry(int sym)
Definition: hypre.cpp:3165
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: hypre.hpp:447
virtual ~HypreParaSails()
Definition: hypre.cpp:3170
The ParaSails preconditioner in hypre.
Definition: hypre.hpp:950
HYPRE_Int GetGlobalNumCols() const
Return the global number of columns.
Definition: hypre.hpp:419
void SetMaxIter(int max_iter)
Definition: hypre.cpp:4237
double * l1_norms
l1 norms of the rows of A
Definition: hypre.hpp:626
Data type sparse matrix.
Definition: sparsemat.hpp:46
void SetOwnerFlags(signed char diag, signed char offd, signed char colmap)
Explicitly set the three ownership flags, see docs for diagOwner etc.
Definition: hypre.hpp:342
void SetLogging(int logging)
Definition: hypre.cpp:2629
Jacobi preconditioner in hypre.
Definition: hypre.hpp:931
virtual ~HypreSolver()
Definition: hypre.cpp:2574
virtual HYPRE_PtrToParSolverFcn SolveFcn() const =0
hypre&#39;s internal Solve function
void SetRelTol(double rel_tol)
Definition: hypre.cpp:4543
void SetKDim(int dim)
Definition: hypre.cpp:2805
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0)
Definition: hypre.cpp:2649
double b
Definition: lissajous.cpp:42
int to_int(const std::string &str)
Convert a string to an int.
Definition: text.hpp:71
Abort on hypre errors (default in base class)
Definition: hypre.hpp:708
void SetLogging(int logging)
Definition: hypre.cpp:2970
void SetPrintLevel(int print_level)
Definition: hypre.hpp:1119
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
FGMRES Setup function.
Definition: hypre.hpp:903
void GetBlocks(Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
Definition: hypre.cpp:929
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
Definition: hypre.cpp:2506
HypreParaSails(MPI_Comm comm)
Definition: hypre.cpp:3079
void SetTol(double tol)
Definition: hypre.cpp:2955
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
PCG Solve function.
Definition: hypre.hpp:806
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:3808
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.cpp:4649
HypreParMatrix * DiscreteCurl(ParFiniteElementSpace *face_fespace, ParFiniteElementSpace *edge_fespace)
Compute the discrete curl matrix between the ND1 and RT0 spaces.
HypreParMatrix * Transpose() const
Returns the transpose of *this.
Definition: hypre.cpp:951
HYPRE_Int Randomize(HYPRE_Int seed)
Set random values.
Definition: hypre.cpp:166
Ignore hypre errors (see e.g. HypreADS)
Definition: hypre.hpp:706
double relax_weight
Damping coefficient (usually &lt;= 1)
Definition: hypre.hpp:610
Parallel smoothers in hypre.
Definition: hypre.hpp:592
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:3561
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:887
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1125
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
ILU Solve function.
Definition: hypre.hpp:1073
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:922
void SetZeroInitialIterate()
non-hypre setting
Definition: hypre.hpp:847
void SetLevelOfFill(HYPRE_Int lev_fill)
Set the fill level for ILU(k); the default is k=1.
Definition: hypre.cpp:3289
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:596
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:942
Dynamic 2D array using row-major layout.
Definition: array.hpp:333
virtual void SetOperator(const Operator &op)
Definition: hypre.cpp:2237
virtual HYPRE_PtrToParSolverFcn SetupFcn() const =0
hypre&#39;s internal Setup function
void SetLogging(int logging)
Definition: hypre.cpp:2810
virtual ~HypreADS()
Definition: hypre.cpp:4094
HypreParMatrix()
An empty matrix to be used as a reference to an existing matrix.
Definition: hypre.cpp:246
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:940
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2624
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2960
HypreParMatrix * GetData()
Definition: hypre.hpp:945
HYPRE_Int * RowPart()
Returns the row partitioning.
Definition: hypre.hpp:366
void SetPolyOptions(int poly_order, double poly_fraction, int eig_est_cg_iter=10)
Set parameters for polynomial smoothing.
Definition: hypre.cpp:2199
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
Definition: hypre.cpp:2193
void BooleanMult(int alpha, const int *x, int beta, int *y)
Definition: hypre.hpp:458
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:70
HypreParMatrix * EliminateCols(const Array< int > &cols)
Definition: hypre.cpp:1391
HypreParVector * X
Definition: hypre.hpp:716
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1211
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:237
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:4323
PCG solver in hypre.
Definition: hypre.hpp:759
GMRES solver in hypre.
Definition: hypre.hpp:817
signed char OwnsOffd() const
Get offd ownership flag.
Definition: hypre.hpp:348
const HYPRE_Int * RowPart() const
Returns the row partitioning (const version)
Definition: hypre.hpp:374
signed char OwnsColMap() const
Get colmap ownership flag.
Definition: hypre.hpp:350
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3839
HypreParVector * X
Definition: hypre.hpp:598
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:4596
void SetNumModes(int num_eigs)
Definition: hypre.cpp:4529
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:977
HypreParMatrix(hypre_ParCSRMatrix *a, bool owner=true)
Converts hypre&#39;s format to HypreParMatrix.
Definition: hypre.hpp:233
virtual ~HypreDiagScale()
Definition: hypre.hpp:946
int GetNumRows() const
Returns the number of rows in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:401
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.hpp:737
HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A + B.
Definition: hypre.cpp:1590
double max_eig_est
Maximal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:632
A class to initialize the size of a Tensor.
Definition: dtensor.hpp:54
void SetRelTol(double rel_tol)
Definition: hypre.cpp:4227
HypreAME(MPI_Comm comm)
Definition: hypre.cpp:4487
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:975
virtual ~HypreIdentity()
Definition: hypre.hpp:927
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
void SetMaxIter(int max_iter)
Definition: hypre.cpp:4553
void SetFIRCoefficients(double max_eig)
Compute window and Chebyshev coefficients for given polynomial order.
Definition: hypre.cpp:2329
void SetRandomSeed(int s)
Definition: hypre.hpp:1337
double a
Definition: lissajous.cpp:41
void EliminateRows(const Array< int > &rows)
Eliminate rows from the diagonal and off-diagonal blocks of the matrix.
Definition: hypre.cpp:1403
signed char OwnsDiag() const
Get diag ownership flag.
Definition: hypre.hpp:346
double * fir_coeffs
Combined coefficients for windowing and Chebyshev polynomials.
Definition: hypre.hpp:639
void Threshold(double threshold=0.0)
Remove values smaller in absolute value than some threshold.
Definition: hypre.cpp:1292
double InnerProduct(HypreParVector *x, HypreParVector *y)
Definition: hypre.cpp:194
HYPRE_Int GlobalSize()
Returns the global number of rows.
Definition: hypre.hpp:122
virtual ~HypreSmoother()
Definition: hypre.cpp:2470
void SetPrintLevel(HYPRE_Int print_level)
Set the print level: 0 = none, 1 = setup, 2 = solve, 3 = setup+solve.
Definition: hypre.cpp:3294
HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col)
Creates vector with given global size and parallel partitioning of the rows/columns given by col...
Definition: hypre.cpp:53
void SetSubSpaceProjector(Operator &proj)
Definition: hypre.hpp:1344
virtual MFEM_DEPRECATED N_Vector ToNVector()
(DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_SERIAL. ...
Definition: vector.hpp:401
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:3441
void SetData(double *_data)
Sets the data of the Vector and the hypre_ParVector to _data.
Definition: hypre.cpp:160
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
Definition: hypre.cpp:1213
MFEM_DEPRECATED void SetZeroInintialIterate()
deprecated: use SetZeroInitialIterate()
Definition: hypre.hpp:787
void SetTaubinOptions(double lambda, double mu, int iter)
Set parameters for Taubin&#39;s lambda-mu method.
Definition: hypre.cpp:2207
HypreParVector * B
Right-hand side and solution vectors.
Definition: hypre.hpp:598
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:3217
int dim
Definition: ex24.cpp:53
void SetMassMatrix(HypreParMatrix &M)
Definition: hypre.cpp:4589
HYPRE_Int M() const
Returns the global number of rows.
Definition: hypre.hpp:380
virtual ~HyprePCG()
Definition: hypre.cpp:2738
int GetNumCols() const
Returns the number of columns in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:408
void SetZeroInitialIterate()
non-hypre setting
Definition: hypre.hpp:790
void SetOperator(Operator &A)
Definition: hypre.cpp:4267
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2639
void GetNumIterations(int &num_iterations)
Definition: hypre.hpp:792
void CopyColStarts()
Definition: hypre.cpp:846
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1016
int relax_times
Number of relaxation sweeps.
Definition: hypre.hpp:608
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
Definition: hypre.cpp:1254
void SetPrintLevel(int logging)
Definition: hypre.cpp:4559
double ParNormlp(const Vector &vec, double p, MPI_Comm comm)
Compute the l_p norm of the Vector which is split without overlap across the given communicator...
Definition: hypre.cpp:205
int eig_est_cg_iter
Number of CG iterations to determine eigenvalue estimates.
Definition: hypre.hpp:630
void SetOwnership(int own)
Sets ownership of the internal hypre_ParVector.
Definition: hypre.hpp:134
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:700
Flexible GMRES solver in hypre.
Definition: hypre.hpp:867
ErrorMode error_mode
How to treat hypre errors.
Definition: hypre.hpp:722
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
ILU Setup function.
Definition: hypre.hpp:1069
HypreParVector & operator=(double d)
Set constant values.
Definition: hypre.cpp:141
const double alpha
Definition: ex15.cpp:336
Vector data type.
Definition: vector.hpp:51
void PrintCommPkg(std::ostream &out=mfem::out) const
Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
Definition: hypre.cpp:1455
HypreParMatrix & operator=(double value)
Initialize all entries with value.
Definition: hypre.hpp:473
hypre_ParCSRMatrix * StealData()
Changes the ownership of the the matrix.
Definition: hypre.cpp:795
void SetPrecondUsageMode(int pcg_mode)
Definition: hypre.cpp:4252
void SetSystemsOptions(int dim, bool order_bynodes=false)
Definition: hypre.cpp:3458
virtual ~HypreParMatrix()
Calls hypre&#39;s destroy function.
Definition: hypre.hpp:547
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:713
ID for class HypreParMatrix.
Definition: operator.hpp:241
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2800
void GetOffd(SparseMatrix &offd, HYPRE_Int *&cmap) const
Get the local off-diagonal block. NOTE: &#39;offd&#39; will not own any data.
Definition: hypre.cpp:923
HYPRE_Int * GetColStarts() const
Return the parallel column partitioning array.
Definition: hypre.hpp:430
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:4638
double omega
SOR parameter (usually in (0,2))
Definition: hypre.hpp:612
HYPRE_Int * GetRowStarts() const
Return the parallel row partitioning array.
Definition: hypre.hpp:425
Base class for solvers.
Definition: operator.hpp:634
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1171
int GetOwnership() const
Gets ownership of the internal hypre_ParVector.
Definition: hypre.hpp:137
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
PCG Setup function.
Definition: hypre.hpp:803
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
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1014
void AbsMult(double a, const Vector &x, double b, Vector &y) const
Computes y = a * |A| * x + b * y, using entry-wise absolute values of matrix A.
Definition: hypre.cpp:1051
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1127
~HypreParVector()
Calls hypre&#39;s destroy function.
Definition: hypre.cpp:176
void SetNumModes(int num_eigs)
Definition: hypre.hpp:1335
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2662
HypreParMatrix * HypreParMatrixFromBlocks(Array2D< HypreParMatrix * > &blocks, Array2D< double > *blockCoeff)
Returns a merged hypre matrix constructed from hypre matrix blocks.
Definition: hypre.cpp:1754
HypreParVector * Z
Definition: hypre.hpp:600
void SetSingularProblem()
Set this option when solving a curl-curl problem with zero mass term.
Definition: hypre.hpp:1166
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:4380
void Read(MPI_Comm comm, const char *fname)
Reads the matrix from a file.
Definition: hypre.cpp:1419
virtual ~HypreFGMRES()
Definition: hypre.cpp:3057
HypreParVector * V
Temporary vectors.
Definition: hypre.hpp:600
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:2187
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
GMRES Setup function.
Definition: hypre.hpp:853
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
int poly_order
Order of the smoothing polynomial.
Definition: hypre.hpp:614
HYPRE_Int * Partitioning()
Returns the parallel row/column partitioning.
Definition: hypre.hpp:119
double lambda
Taubin&#39;s lambda-mu method parameters.
Definition: hypre.hpp:621
MFEM_DEPRECATED void SetZeroInintialIterate()
deprecated: use SetZeroInitialIterate()
Definition: hypre.hpp:894
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1213