MFEM  v4.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, 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_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 #ifdef MFEM_USE_SUNDIALS
38 #include <nvector/nvector_parhyp.h>
39 #endif
40 
41 namespace mfem
42 {
43 
44 class ParFiniteElementSpace;
45 class HypreParMatrix;
46 
47 namespace internal
48 {
49 
50 template <typename int_type>
51 inline int to_int(int_type i)
52 {
53  MFEM_ASSERT(int_type(int(i)) == i, "overflow converting int_type to int");
54  return int(i);
55 }
56 
57 // Specialization for to_int(int)
58 template <> inline int to_int(int i) { return i; }
59 
60 // Convert a HYPRE_Int to int
61 #ifdef HYPRE_BIGINT
62 template <>
63 inline int to_int(HYPRE_Int i)
64 {
65  MFEM_ASSERT(HYPRE_Int(int(i)) == i, "overflow converting HYPRE_Int to int");
66  return int(i);
67 }
68 #endif
69 
70 }
71 
72 /// Wrapper for hypre's parallel vector class
73 class HypreParVector : public Vector
74 {
75 private:
76  int own_ParVector;
77 
78  /// The actual object
79  hypre_ParVector *x;
80 
81  friend class HypreParMatrix;
82 
83  // Set Vector::data and Vector::size from *x
84  inline void _SetDataAndSize_();
85 
86 public:
87  /** @brief Creates vector with given global size and parallel partitioning of
88  the rows/columns given by @a col. */
89  /** @anchor hypre_partitioning_descr
90  The partitioning is defined in one of two ways depending on the
91  configuration of HYPRE:
92  1. If HYPRE_AssumedPartitionCheck() returns true (the default),
93  then col is of length 2 and the local processor owns columns
94  [col[0],col[1]).
95  2. If HYPRE_AssumedPartitionCheck() returns false, then col is of
96  length (number of processors + 1) and processor P owns columns
97  [col[P],col[P+1]) i.e. each processor has a copy of the same col
98  array. */
99  HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col);
100  /** @brief Creates vector with given global size, partitioning of the
101  columns, and data. */
102  /** The data must be allocated and destroyed outside. If @a _data is NULL, a
103  dummy vector without a valid data array will be created. See @ref
104  hypre_partitioning_descr "here" for a description of the @a col array. */
105  HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, double *_data,
106  HYPRE_Int *col);
107  /// Creates vector compatible with y
108  HypreParVector(const HypreParVector &y);
109  /// Creates vector compatible with (i.e. in the domain of) A or A^T
110  explicit HypreParVector(const HypreParMatrix &A, int transpose = 0);
111  /// Creates vector wrapping y
112  explicit HypreParVector(HYPRE_ParVector y);
113  /// Create a true dof parallel vector on a given ParFiniteElementSpace
114  explicit HypreParVector(ParFiniteElementSpace *pfes);
115 
116  /// MPI communicator
117  MPI_Comm GetComm() { return x->comm; }
118 
119  /// Returns the parallel row/column partitioning
120  /** See @ref hypre_partitioning_descr "here" for a description of the
121  partitioning array. */
122  inline HYPRE_Int *Partitioning() { return x->partitioning; }
123 
124  /// Returns the global number of rows
125  inline HYPRE_Int GlobalSize() { return x->global_size; }
126 
127  /// Typecasting to hypre's hypre_ParVector*
128  operator hypre_ParVector*() const { return x; }
129 #ifndef HYPRE_PAR_VECTOR_STRUCT
130  /// Typecasting to hypre's HYPRE_ParVector, a.k.a. void *
131  operator HYPRE_ParVector() const { return (HYPRE_ParVector) x; }
132 #endif
133  /// Changes the ownership of the the vector
134  hypre_ParVector *StealParVector() { own_ParVector = 0; return x; }
135 
136  /// Sets ownership of the internal hypre_ParVector
137  void SetOwnership(int own) { own_ParVector = own; }
138 
139  /// Gets ownership of the internal hypre_ParVector
140  int GetOwnership() const { return own_ParVector; }
141 
142  /// Returns the global vector in each processor
143  Vector* GlobalVector() const;
144 
145  /// Set constant values
146  HypreParVector& operator= (double d);
147  /// Define '=' for hypre vectors.
149 
150  /// Sets the data of the Vector and the hypre_ParVector to @a _data.
151  /** Must be used only for HypreParVector%s that do not own the data,
152  e.g. created with the constructor:
153  HypreParVector(MPI_Comm, HYPRE_Int, double *, HYPRE_Int *). */
154  void SetData(double *_data);
155 
156  /// Set random values
157  HYPRE_Int Randomize(HYPRE_Int seed);
158 
159  /// Prints the locally owned rows in parallel
160  void Print(const char *fname) const;
161 
162  /// Calls hypre's destroy function
163  ~HypreParVector();
164 
165 #ifdef MFEM_USE_SUNDIALS
166  /// Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_PARHYP.
167  /** The returned N_Vector must be destroyed by the caller. */
168  virtual N_Vector ToNVector() { return N_VMake_ParHyp(x); }
169 
170  /** @brief Update an existing wrapper SUNDIALS N_Vector of type
171  SUNDIALS_NVEC_PARHYP to point to this Vector. */
172  virtual void ToNVector(N_Vector &nv);
173 #endif
174 };
175 
176 /// Returns the inner product of x and y
177 double InnerProduct(HypreParVector &x, HypreParVector &y);
178 double InnerProduct(HypreParVector *x, HypreParVector *y);
179 
180 
181 /** @brief Compute the l_p norm of the Vector which is split without overlap
182  across the given communicator. */
183 double ParNormlp(const Vector &vec, double p, MPI_Comm comm);
184 
185 
186 /// Wrapper for hypre's ParCSR matrix class
187 class HypreParMatrix : public Operator
188 {
189 private:
190  /// The actual object
191  hypre_ParCSRMatrix *A;
192 
193  /// Auxiliary vectors for typecasting
194  mutable HypreParVector *X, *Y;
195 
196  // Flags indicating ownership of A->diag->{i,j,data}, A->offd->{i,j,data},
197  // and A->col_map_offd.
198  // The possible values for diagOwner are:
199  // -1: no special treatment of A->diag (default)
200  // 0: prevent hypre from destroying A->diag->{i,j,data}
201  // 1: same as 0, plus take ownership of A->diag->{i,j}
202  // 2: same as 0, plus take ownership of A->diag->data
203  // 3: same as 0, plus take ownership of A->diag->{i,j,data}
204  // The same values and rules apply to offdOwner and A->offd.
205  // The possible values for colMapOwner are:
206  // -1: no special treatment of A->col_map_offd (default)
207  // 0: prevent hypre from destroying A->col_map_offd
208  // 1: same as 0, plus take ownership of A->col_map_offd
209  // All owned arrays are destroyed with 'delete []'.
210  signed char diagOwner, offdOwner, colMapOwner;
211 
212  // Does the object own the pointer A?
213  signed char ParCSROwner;
214 
215  // Initialize with defaults. Does not initialize inherited members.
216  void Init();
217 
218  // Delete all owned data. Does not perform re-initialization with defaults.
219  void Destroy();
220 
221  // Copy (shallow/deep, based on HYPRE_BIGINT) the I and J arrays from csr to
222  // hypre_csr. Shallow copy the data. Return the appropriate ownership flag.
223  static char CopyCSR(SparseMatrix *csr, hypre_CSRMatrix *hypre_csr);
224  // Copy (shallow or deep, based on HYPRE_BIGINT) the I and J arrays from
225  // bool_csr to hypre_csr. Allocate the data array and set it to all ones.
226  // Return the appropriate ownership flag.
227  static char CopyBoolCSR(Table *bool_csr, hypre_CSRMatrix *hypre_csr);
228 
229  // Copy the j array of a hypre_CSRMatrix to the given J array, converting
230  // the indices from HYPRE_Int to int.
231  static void CopyCSR_J(hypre_CSRMatrix *hypre_csr, int *J);
232 
233 public:
234  /// An empty matrix to be used as a reference to an existing matrix
235  HypreParMatrix();
236 
237  /// Converts hypre's format to HypreParMatrix
238  /** If @a owner is false, ownership of @a a is not transferred */
239  explicit HypreParMatrix(hypre_ParCSRMatrix *a, bool owner = true)
240  {
241  Init();
242  A = a;
243  if (!owner) { ParCSROwner = 0; }
244  height = GetNumRows();
245  width = GetNumCols();
246  }
247 
248  /// Creates block-diagonal square parallel matrix.
249  /** Diagonal is given by @a diag which must be in CSR format (finalized). The
250  new HypreParMatrix does not take ownership of any of the input arrays.
251  See @ref hypre_partitioning_descr "here" for a description of the row
252  partitioning array @a row_starts.
253 
254  @warning The ordering of the columns in each row in @a *diag may be
255  changed by this constructor to ensure that the first entry in each row is
256  the diagonal one. This is expected by most hypre functions. */
257  HypreParMatrix(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *row_starts,
258  SparseMatrix *diag); // constructor with 4 arguments, v1
259 
260  /// Creates block-diagonal rectangular parallel matrix.
261  /** Diagonal is given by @a diag which must be in CSR format (finalized). The
262  new HypreParMatrix does not take ownership of any of the input arrays.
263  See @ref hypre_partitioning_descr "here" for a description of the
264  partitioning arrays @a row_starts and @a col_starts. */
265  HypreParMatrix(MPI_Comm comm, HYPRE_Int global_num_rows,
266  HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
267  HYPRE_Int *col_starts,
268  SparseMatrix *diag); // constructor with 6 arguments, v1
269 
270  /// Creates general (rectangular) parallel matrix.
271  /** The new HypreParMatrix does not take ownership of any of the input
272  arrays. See @ref hypre_partitioning_descr "here" for a description of the
273  partitioning arrays @a row_starts and @a col_starts. */
274  HypreParMatrix(MPI_Comm comm, HYPRE_Int global_num_rows,
275  HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
276  HYPRE_Int *col_starts, SparseMatrix *diag, SparseMatrix *offd,
277  HYPRE_Int *cmap); // constructor with 8 arguments
278 
279  /// Creates general (rectangular) parallel matrix.
280  /** The new HypreParMatrix takes ownership of all input arrays, except
281  @a col_starts and @a row_starts. See @ref hypre_partitioning_descr "here"
282  for a description of the partitioning arrays @a row_starts and @a
283  col_starts. */
284  HypreParMatrix(MPI_Comm comm,
285  HYPRE_Int global_num_rows, HYPRE_Int global_num_cols,
286  HYPRE_Int *row_starts, HYPRE_Int *col_starts,
287  HYPRE_Int *diag_i, HYPRE_Int *diag_j, double *diag_data,
288  HYPRE_Int *offd_i, HYPRE_Int *offd_j, double *offd_data,
289  HYPRE_Int offd_num_cols,
290  HYPRE_Int *offd_col_map); // constructor with 13 arguments
291 
292  /// Creates a parallel matrix from SparseMatrix on processor 0.
293  /** See @ref hypre_partitioning_descr "here" for a description of the
294  partitioning arrays @a row_starts and @a col_starts. */
295  HypreParMatrix(MPI_Comm comm, HYPRE_Int *row_starts, HYPRE_Int *col_starts,
296  SparseMatrix *a); // constructor with 4 arguments, v2
297 
298  /// Creates boolean block-diagonal rectangular parallel matrix.
299  /** The new HypreParMatrix does not take ownership of any of the input
300  arrays. See @ref hypre_partitioning_descr "here" for a description of the
301  partitioning arrays @a row_starts and @a col_starts. */
302  HypreParMatrix(MPI_Comm comm, HYPRE_Int global_num_rows,
303  HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
304  HYPRE_Int *col_starts,
305  Table *diag); // constructor with 6 arguments, v2
306 
307  /// Creates boolean rectangular parallel matrix.
308  /** The new HypreParMatrix takes ownership of the arrays @a i_diag,
309  @a j_diag, @a i_offd, @a j_offd, and @a cmap; does not take ownership of
310  the arrays @a row and @a col. See @ref hypre_partitioning_descr "here"
311  for a description of the partitioning arrays @a row and @a col. */
312  HypreParMatrix(MPI_Comm comm, int id, int np, HYPRE_Int *row, HYPRE_Int *col,
313  HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd,
314  HYPRE_Int *j_offd, HYPRE_Int *cmap,
315  HYPRE_Int cmap_size); // constructor with 11 arguments
316 
317  /** @brief Creates a general parallel matrix from a local CSR matrix on each
318  processor described by the @a I, @a J and @a data arrays. */
319  /** The local matrix should be of size (local) @a nrows by (global)
320  @a glob_ncols. The new parallel matrix contains copies of all input
321  arrays (so they can be deleted). See @ref hypre_partitioning_descr "here"
322  for a description of the partitioning arrays @a rows and @a cols. */
323  HypreParMatrix(MPI_Comm comm, int nrows, HYPRE_Int glob_nrows,
324  HYPRE_Int glob_ncols, int *I, HYPRE_Int *J,
325  double *data, HYPRE_Int *rows,
326  HYPRE_Int *cols); // constructor with 9 arguments
327 
328  /// Make this HypreParMatrix a reference to 'master'
329  void MakeRef(const HypreParMatrix &master);
330 
331  /// MPI communicator
332  MPI_Comm GetComm() const { return A->comm; }
333 
334  /// Typecasting to hypre's hypre_ParCSRMatrix*
335  operator hypre_ParCSRMatrix*() const { return A; }
336 #ifndef HYPRE_PAR_CSR_MATRIX_STRUCT
337  /// Typecasting to hypre's HYPRE_ParCSRMatrix, a.k.a. void *
338  operator HYPRE_ParCSRMatrix() { return (HYPRE_ParCSRMatrix) A; }
339 #endif
340  /// Changes the ownership of the the matrix
341  hypre_ParCSRMatrix* StealData();
342 
343  /// Explicitly set the three ownership flags, see docs for diagOwner etc.
344  void SetOwnerFlags(signed char diag, signed char offd, signed char colmap)
345  { diagOwner = diag, offdOwner = offd, colMapOwner = colmap; }
346 
347  /// Get diag ownership flag
348  signed char OwnsDiag() const { return diagOwner; }
349  /// Get offd ownership flag
350  signed char OwnsOffd() const { return offdOwner; }
351  /// Get colmap ownership flag
352  signed char OwnsColMap() const { return colMapOwner; }
353 
354  /** If the HypreParMatrix does not own the row-starts array, make a copy of
355  it that the HypreParMatrix will own. If the col-starts array is the same
356  as the row-starts array, col-starts is also replaced. */
357  void CopyRowStarts();
358  /** If the HypreParMatrix does not own the col-starts array, make a copy of
359  it that the HypreParMatrix will own. If the row-starts array is the same
360  as the col-starts array, row-starts is also replaced. */
361  void CopyColStarts();
362 
363  /// Returns the global number of nonzeros
364  inline HYPRE_Int NNZ() const { return A->num_nonzeros; }
365  /// Returns the row partitioning
366  /** See @ref hypre_partitioning_descr "here" for a description of the
367  partitioning array. */
368  inline HYPRE_Int *RowPart() { return A->row_starts; }
369  /// Returns the column partitioning
370  /** See @ref hypre_partitioning_descr "here" for a description of the
371  partitioning array. */
372  inline HYPRE_Int *ColPart() { return A->col_starts; }
373  /// Returns the row partitioning (const version)
374  /** See @ref hypre_partitioning_descr "here" for a description of the
375  partitioning array. */
376  inline const HYPRE_Int *RowPart() const { return A->row_starts; }
377  /// Returns the column partitioning (const version)
378  /** See @ref hypre_partitioning_descr "here" for a description of the
379  partitioning array. */
380  inline const HYPRE_Int *ColPart() const { return A->col_starts; }
381  /// Returns the global number of rows
382  inline HYPRE_Int M() const { return A->global_num_rows; }
383  /// Returns the global number of columns
384  inline HYPRE_Int N() const { return A->global_num_cols; }
385 
386  /// Get the local diagonal of the matrix.
387  void GetDiag(Vector &diag) const;
388  /// Get the local diagonal block. NOTE: 'diag' will not own any data.
389  void GetDiag(SparseMatrix &diag) const;
390  /// Get the local off-diagonal block. NOTE: 'offd' will not own any data.
391  void GetOffd(SparseMatrix &offd, HYPRE_Int* &cmap) const;
392 
393  /** Split the matrix into M x N equally sized blocks of parallel matrices.
394  The size of 'blocks' must already be set to M x N. */
395  void GetBlocks(Array2D<HypreParMatrix*> &blocks,
396  bool interleaved_rows = false,
397  bool interleaved_cols = false) const;
398 
399  /// Returns the transpose of *this
400  HypreParMatrix * Transpose() const;
401 
402  /// Returns the number of rows in the diagonal block of the ParCSRMatrix
403  int GetNumRows() const
404  {
405  return internal::to_int(
406  hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(A)));
407  }
408 
409  /// Returns the number of columns in the diagonal block of the ParCSRMatrix
410  int GetNumCols() const
411  {
412  return internal::to_int(
413  hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(A)));
414  }
415 
416  /// Return the global number of rows
417  HYPRE_Int GetGlobalNumRows() const
418  { return hypre_ParCSRMatrixGlobalNumRows(A); }
419 
420  /// Return the global number of columns
421  HYPRE_Int GetGlobalNumCols() const
422  { return hypre_ParCSRMatrixGlobalNumCols(A); }
423 
424  /// Return the parallel row partitioning array.
425  /** See @ref hypre_partitioning_descr "here" for a description of the
426  partitioning array. */
427  HYPRE_Int *GetRowStarts() const { return hypre_ParCSRMatrixRowStarts(A); }
428 
429  /// Return the parallel column partitioning array.
430  /** See @ref hypre_partitioning_descr "here" for a description of the
431  partitioning array. */
432  HYPRE_Int *GetColStarts() const { return hypre_ParCSRMatrixColStarts(A); }
433 
434  /// Computes y = alpha * A * x + beta * y
435  HYPRE_Int Mult(HypreParVector &x, HypreParVector &y,
436  double alpha = 1.0, double beta = 0.0);
437  /// Computes y = alpha * A * x + beta * y
438  HYPRE_Int Mult(HYPRE_ParVector x, HYPRE_ParVector y,
439  double alpha = 1.0, double beta = 0.0);
440  /// Computes y = alpha * A^t * x + beta * y
442  double alpha = 1.0, double beta = 0.0);
443 
444  void Mult(double a, const Vector &x, double b, Vector &y) const;
445  void MultTranspose(double a, const Vector &x, double b, Vector &y) const;
446 
447  virtual void Mult(const Vector &x, Vector &y) const
448  { Mult(1.0, x, 0.0, y); }
449  virtual void MultTranspose(const Vector &x, Vector &y) const
450  { MultTranspose(1.0, x, 0.0, y); }
451 
452  /** The "Boolean" analog of y = alpha * A * x + beta * y, where elements in
453  the sparsity pattern of the matrix are treated as "true". */
454  void BooleanMult(int alpha, const int *x, int beta, int *y)
455  {
456  internal::hypre_ParCSRMatrixBooleanMatvec(A, alpha, const_cast<int*>(x),
457  beta, y);
458  }
459 
460  /** The "Boolean" analog of y = alpha * A^T * x + beta * y, where elements in
461  the sparsity pattern of the matrix are treated as "true". */
462  void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
463  {
464  internal::hypre_ParCSRMatrixBooleanMatvecT(A, alpha, const_cast<int*>(x),
465  beta, y);
466  }
467 
468  /// Initialize all entries with value.
469  HypreParMatrix &operator=(double value)
470  { internal::hypre_ParCSRMatrixSetConstantValues(A, value); return *this; }
471 
472  /** Perform the operation `*this += B`, assuming that both matrices use the
473  same row and column partitions and the same col_map_offd arrays, or B has
474  an empty off-diagonal block. We also assume that the sparsity pattern of
475  `*this` contains that of `B`. */
476  HypreParMatrix &operator+=(const HypreParMatrix &B) { return Add(1.0, B); }
477 
478  /** Perform the operation `*this += beta*B`, assuming that both matrices use
479  the same row and column partitions and the same col_map_offd arrays, or
480  B has an empty off-diagonal block. We also assume that the sparsity
481  pattern of `*this` contains that of `B`. For a more general case consider
482  the stand-alone function ParAdd described below. */
483  HypreParMatrix &Add(const double beta, const HypreParMatrix &B)
484  {
485  MFEM_VERIFY(internal::hypre_ParCSRMatrixSum(A, beta, B.A) == 0,
486  "error in hypre_ParCSRMatrixSum");
487  return *this;
488  }
489 
490  /** @brief Multiply the HypreParMatrix on the left by a block-diagonal
491  parallel matrix @a D and return the result as a new HypreParMatrix. */
492  /** If @a D has a different number of rows than @a A (this matrix), @a D's
493  row starts array needs to be given (as returned by the methods
494  GetDofOffsets/GetTrueDofOffsets of ParFiniteElementSpace). The new
495  matrix @a D*A uses copies of the row- and column-starts arrays, so "this"
496  matrix and @a row_starts can be deleted.
497  @note This operation is local and does not require communication. */
499  HYPRE_Int* row_starts = NULL) const;
500 
501  /// Scale the local row i by s(i).
502  void ScaleRows(const Vector & s);
503  /// Scale the local row i by 1./s(i)
504  void InvScaleRows(const Vector & s);
505  /// Scale all entries by s: A_scaled = s*A.
506  void operator*=(double s);
507 
508  /// Remove values smaller in absolute value than some threshold
509  void Threshold(double threshold = 0.0);
510 
511  /// If a row contains only zeros, set its diagonal to 1.
512  void EliminateZeroRows() { hypre_ParCSRMatrixFixZeroRows(A); }
513 
514  /** Eliminate rows and columns from the matrix, and rows from the vector B.
515  Modify B with the BC values in X. */
516  void EliminateRowsCols(const Array<int> &rows_cols, const HypreParVector &X,
517  HypreParVector &B);
518 
519  /** Eliminate rows and columns from the matrix and store the eliminated
520  elements in a new matrix Ae (returned), so that the modified matrix and
521  Ae sum to the original matrix. */
522  HypreParMatrix* EliminateRowsCols(const Array<int> &rows_cols);
523 
524  /// Prints the locally owned rows in parallel
525  void Print(const char *fname, HYPRE_Int offi = 0, HYPRE_Int offj = 0);
526  /// Reads the matrix from a file
527  void Read(MPI_Comm comm, const char *fname);
528  /// Read a matrix saved as a HYPRE_IJMatrix
529  void Read_IJMatrix(MPI_Comm comm, const char *fname);
530 
531  /// Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
532  void PrintCommPkg(std::ostream &out = mfem::out) const;
533 
534  /// Calls hypre's destroy function
535  virtual ~HypreParMatrix() { Destroy(); }
536 
537  Type GetType() const { return Hypre_ParCSR; }
538 };
539 
540 /** @brief Return a new matrix `C = alpha*A + beta*B`, assuming that both `A`
541  and `B` use the same row and column partitions and the same `col_map_offd`
542  arrays. */
543 HypreParMatrix *Add(double alpha, const HypreParMatrix &A,
544  double beta, const HypreParMatrix &B);
545 
546 /// Returns the matrix A * B
547 HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B);
548 /// Returns the matrix A + B
549 /** It is assumed that both matrices use the same row and column partitions and
550  the same col_map_offd arrays. */
551 HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B);
552 
553 /// Returns the matrix P^t * A * P
554 HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P);
555 /// Returns the matrix Rt^t * A * P
556 HypreParMatrix * RAP(const HypreParMatrix * Rt, const HypreParMatrix *A,
557  const HypreParMatrix *P);
558 
559 /** Eliminate essential BC specified by 'ess_dof_list' from the solution X to
560  the r.h.s. B. Here A is a matrix with eliminated BC, while Ae is such that
561  (A+Ae) is the original (Neumann) matrix before elimination. */
562 void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae,
563  const Array<int> &ess_dof_list, const Vector &X, Vector &B);
564 
565 
566 /// Parallel smoothers in hypre
567 class HypreSmoother : public Solver
568 {
569 protected:
570  /// The linear system matrix
572  /// Right-hand side and solution vectors
573  mutable HypreParVector *B, *X;
574  /// Temporary vectors
575  mutable HypreParVector *V, *Z;
576  /// FIR Filter Temporary Vectors
577  mutable HypreParVector *X0, *X1;
578 
579  /** Smoother type from hypre_ParCSRRelax() in ams.c plus extensions, see the
580  enumeration Type below. */
581  int type;
582  /// Number of relaxation sweeps
584  /// Damping coefficient (usually <= 1)
585  double relax_weight;
586  /// SOR parameter (usually in (0,2))
587  double omega;
588  /// Order of the smoothing polynomial
590  /// Fraction of spectrum to smooth for polynomial relaxation
592  /// Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}
594 
595  /// Taubin's lambda-mu method parameters
596  double lambda;
597  double mu;
599 
600  /// l1 norms of the rows of A
601  double *l1_norms;
602  /// If set, take absolute values of the computed l1_norms
604  /// Maximal eigenvalue estimate for polynomial smoothing
605  double max_eig_est;
606  /// Minimal eigenvalue estimate for polynomial smoothing
607  double min_eig_est;
608  /// Parameters for windowing function of FIR filter
609  double window_params[3];
610 
611  /// Combined coefficients for windowing and Chebyshev polynomials.
612  double* fir_coeffs;
613 
614 public:
615  /** Hypre smoother types:
616  0 = Jacobi
617  1 = l1-scaled Jacobi
618  2 = l1-scaled block Gauss-Seidel/SSOR
619  4 = truncated l1-scaled block Gauss-Seidel/SSOR
620  5 = lumped Jacobi
621  6 = Gauss-Seidel
622  16 = Chebyshev
623  1001 = Taubin polynomial smoother
624  1002 = FIR polynomial smoother. */
625  enum Type { Jacobi = 0, l1Jacobi = 1, l1GS = 2, l1GStr = 4, lumpedJacobi = 5,
626  GS = 6, Chebyshev = 16, Taubin = 1001, FIR = 1002
627  };
628 
629  HypreSmoother();
630 
632  int relax_times = 1, double relax_weight = 1.0,
633  double omega = 1.0, int poly_order = 2,
634  double poly_fraction = .3);
635 
636  /// Set the relaxation type and number of sweeps
638  /// Set SOR-related parameters
639  void SetSOROptions(double relax_weight, double omega);
640  /// Set parameters for polynomial smoothing
641  void SetPolyOptions(int poly_order, double poly_fraction);
642  /// Set parameters for Taubin's lambda-mu method
643  void SetTaubinOptions(double lambda, double mu, int iter);
644 
645  /// Convenience function for setting canonical windowing parameters
646  void SetWindowByName(const char* window_name);
647  /// Set parameters for windowing function for FIR smoother.
648  void SetWindowParameters(double a, double b, double c);
649  /// Compute window and Chebyshev coefficients for given polynomial order.
650  void SetFIRCoefficients(double max_eig);
651 
652  /// After computing l1-norms, replace them with their absolute values.
653  /** By default, the l1-norms take their sign from the corresponding diagonal
654  entries in the associated matrix. */
655  void SetPositiveDiagonal(bool pos = true) { pos_l1_norms = pos; }
656 
657  /** Set/update the associated operator. Must be called after setting the
658  HypreSmoother type and options. */
659  virtual void SetOperator(const Operator &op);
660 
661  /// Relax the linear system Ax=b
662  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
663  virtual void Mult(const Vector &b, Vector &x) const;
664 
665  virtual ~HypreSmoother();
666 };
667 
668 
669 /// Abstract class for hypre's solvers and preconditioners
670 class HypreSolver : public Solver
671 {
672 protected:
673  /// The linear system matrix
675 
676  /// Right-hand side and solution vector
677  mutable HypreParVector *B, *X;
678 
679  /// Was hypre's Setup function called already?
680  mutable int setup_called;
681 
682 public:
683  HypreSolver();
684 
686 
687  /// Typecast to HYPRE_Solver -- return the solver
688  virtual operator HYPRE_Solver() const = 0;
689 
690  /// hypre's internal Setup function
691  virtual HYPRE_PtrToParSolverFcn SetupFcn() const = 0;
692  /// hypre's internal Solve function
693  virtual HYPRE_PtrToParSolverFcn SolveFcn() const = 0;
694 
695  virtual void SetOperator(const Operator &op)
696  { mfem_error("HypreSolvers do not support SetOperator!"); }
697 
698  /// Solve the linear system Ax=b
699  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
700  virtual void Mult(const Vector &b, Vector &x) const;
701 
702  virtual ~HypreSolver();
703 };
704 
705 /// PCG solver in hypre
706 class HyprePCG : public HypreSolver
707 {
708 private:
709  HYPRE_Solver pcg_solver;
710 
711 public:
713 
714  void SetTol(double tol);
715  void SetMaxIter(int max_iter);
716  void SetLogging(int logging);
717  void SetPrintLevel(int print_lvl);
718 
719  /// Set the hypre solver to be used as a preconditioner
720  void SetPreconditioner(HypreSolver &precond);
721 
722  /** Use the L2 norm of the residual for measuring PCG convergence, plus
723  (optionally) 1) periodically recompute true residuals from scratch; and
724  2) enable residual-based stopping criteria. */
725  void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0);
726 
727  /// non-hypre setting
729 
730  void GetNumIterations(int &num_iterations)
731  {
732  HYPRE_Int num_it;
733  HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_it);
734  num_iterations = internal::to_int(num_it);
735  }
736 
737  /// The typecast to HYPRE_Solver returns the internal pcg_solver
738  virtual operator HYPRE_Solver() const { return pcg_solver; }
739 
740  /// PCG Setup function
741  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
742  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSetup; }
743  /// PCG Solve function
744  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
745  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSolve; }
746 
747  /// Solve Ax=b with hypre's PCG
748  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
749  using HypreSolver::Mult;
750 
751  virtual ~HyprePCG();
752 };
753 
754 /// GMRES solver in hypre
755 class HypreGMRES : public HypreSolver
756 {
757 private:
758  HYPRE_Solver gmres_solver;
759 
760 public:
762 
763  void SetTol(double tol);
764  void SetMaxIter(int max_iter);
765  void SetKDim(int dim);
766  void SetLogging(int logging);
767  void SetPrintLevel(int print_lvl);
768 
769  /// Set the hypre solver to be used as a preconditioner
770  void SetPreconditioner(HypreSolver &precond);
771 
772  /// non-hypre setting
774 
775  /// The typecast to HYPRE_Solver returns the internal gmres_solver
776  virtual operator HYPRE_Solver() const { return gmres_solver; }
777 
778  /// GMRES Setup function
779  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
780  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSetup; }
781  /// GMRES Solve function
782  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
783  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSolve; }
784 
785  /// Solve Ax=b with hypre's GMRES
786  virtual void Mult (const HypreParVector &b, HypreParVector &x) const;
787  using HypreSolver::Mult;
788 
789  virtual ~HypreGMRES();
790 };
791 
792 /// The identity operator as a hypre solver
794 {
795 public:
796  virtual operator HYPRE_Solver() const { return NULL; }
797 
798  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
799  { return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentitySetup; }
800  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
801  { return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentity; }
802 
803  virtual ~HypreIdentity() { }
804 };
805 
806 /// Jacobi preconditioner in hypre
808 {
809 public:
812  virtual operator HYPRE_Solver() const { return NULL; }
813 
814  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
815  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScaleSetup; }
816  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
817  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScale; }
818 
819  HypreParMatrix* GetData() { return A; }
820  virtual ~HypreDiagScale() { }
821 };
822 
823 /// The ParaSails preconditioner in hypre
825 {
826 private:
827  HYPRE_Solver sai_precond;
828 
829 public:
831 
832  void SetSymmetry(int sym);
833 
834  /// The typecast to HYPRE_Solver returns the internal sai_precond
835  virtual operator HYPRE_Solver() const { return sai_precond; }
836 
837  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
838  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSetup; }
839  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
840  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSolve; }
841 
842  virtual ~HypreParaSails();
843 };
844 
845 /** The Euclid preconditioner in Hypre
846 
847  Euclid implements the Parallel Incomplete LU factorization technique. For
848  more information see:
849 
850  "A Scalable Parallel Algorithm for Incomplete Factor Preconditioning" by
851  David Hysom and Alex Pothen, https://doi.org/10.1137/S1064827500376193
852 */
853 class HypreEuclid : public HypreSolver
854 {
855 private:
856  HYPRE_Solver euc_precond;
857 
858 public:
860 
861  /// The typecast to HYPRE_Solver returns the internal euc_precond
862  virtual operator HYPRE_Solver() const { return euc_precond; }
863 
864  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
865  { return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSetup; }
866  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
867  { return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSolve; }
868 
869  virtual ~HypreEuclid();
870 };
871 
872 /// The BoomerAMG solver in hypre
874 {
875 private:
876  HYPRE_Solver amg_precond;
877 
878  /// Rigid body modes
880 
881  /// Finite element space for elasticity problems, see SetElasticityOptions()
882  ParFiniteElementSpace *fespace;
883 
884  /// Recompute the rigid-body modes vectors (in the rbms array)
885  void RecomputeRBMs();
886 
887  /// Default, generally robust, BoomerAMG options
888  void SetDefaultOptions();
889 
890  // If amg_precond is NULL, allocates it and sets default options.
891  // Otherwise saves the options from amg_precond, destroys it, allocates a new
892  // one, and sets its options to the saved values.
893  void ResetAMGPrecond();
894 
895 public:
896  HypreBoomerAMG();
897 
899 
900  virtual void SetOperator(const Operator &op);
901 
902  /** More robust options for systems, such as elasticity. Note that BoomerAMG
903  assumes Ordering::byVDIM in the finite element space used to generate the
904  matrix A. */
905  void SetSystemsOptions(int dim);
906 
907  /** A special elasticity version of BoomerAMG that takes advantage of
908  geometric rigid body modes and could perform better on some problems, see
909  "Improving algebraic multigrid interpolation operators for linear
910  elasticity problems", Baker, Kolev, Yang, NLAA 2009, DOI:10.1002/nla.688.
911  As with SetSystemsOptions(), this solver assumes Ordering::byVDIM. */
913 
914  void SetPrintLevel(int print_level)
915  { HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level); }
916 
917  /// The typecast to HYPRE_Solver returns the internal amg_precond
918  virtual operator HYPRE_Solver() const { return amg_precond; }
919 
920  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
921  { return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSetup; }
922  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
923  { return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSolve; }
924 
925  virtual ~HypreBoomerAMG();
926 };
927 
928 /// Compute the discrete gradient matrix between the nodal linear and ND1 spaces
929 HypreParMatrix* DiscreteGrad(ParFiniteElementSpace *edge_fespace,
930  ParFiniteElementSpace *vert_fespace);
931 /// Compute the discrete curl matrix between the ND1 and RT0 spaces
932 HypreParMatrix* DiscreteCurl(ParFiniteElementSpace *face_fespace,
933  ParFiniteElementSpace *edge_fespace);
934 
935 /// The Auxiliary-space Maxwell Solver in hypre
936 class HypreAMS : public HypreSolver
937 {
938 private:
939  HYPRE_Solver ams;
940 
941  /// Vertex coordinates
942  HypreParVector *x, *y, *z;
943  /// Discrete gradient matrix
944  HypreParMatrix *G;
945  /// Nedelec interpolation matrix and its components
946  HypreParMatrix *Pi, *Pix, *Piy, *Piz;
947 
948 public:
950 
951  void SetPrintLevel(int print_lvl);
952 
953  /// Set this option when solving a curl-curl problem with zero mass term
954  void SetSingularProblem() { HYPRE_AMSSetBetaPoissonMatrix(ams, NULL); }
955 
956  /// The typecast to HYPRE_Solver returns the internal ams object
957  virtual operator HYPRE_Solver() const { return ams; }
958 
959  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
960  { return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSetup; }
961  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
962  { return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSolve; }
963 
964  virtual ~HypreAMS();
965 };
966 
967 /// The Auxiliary-space Divergence Solver in hypre
968 class HypreADS : public HypreSolver
969 {
970 private:
971  HYPRE_Solver ads;
972 
973  /// Vertex coordinates
974  HypreParVector *x, *y, *z;
975  /// Discrete gradient matrix
976  HypreParMatrix *G;
977  /// Discrete curl matrix
978  HypreParMatrix *C;
979  /// Nedelec interpolation matrix and its components
980  HypreParMatrix *ND_Pi, *ND_Pix, *ND_Piy, *ND_Piz;
981  /// Raviart-Thomas interpolation matrix and its components
982  HypreParMatrix *RT_Pi, *RT_Pix, *RT_Piy, *RT_Piz;
983 
984 public:
986 
987  void SetPrintLevel(int print_lvl);
988 
989  /// The typecast to HYPRE_Solver returns the internal ads object
990  virtual operator HYPRE_Solver() const { return ads; }
991 
992  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
993  { return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSetup; }
994  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
995  { return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSolve; }
996 
997  virtual ~HypreADS();
998 };
999 
1000 /** LOBPCG eigenvalue solver in hypre
1001 
1002  The Locally Optimal Block Preconditioned Conjugate Gradient (LOBPCG)
1003  eigenvalue solver is designed to find the lowest eigenmodes of the
1004  generalized eigenvalue problem:
1005  A x = lambda M x
1006  where A is symmetric, potentially indefinite and M is symmetric positive
1007  definite. The eigenvectors are M-orthonormal, meaning that
1008  x^T M x = 1 and x^T M y = 0,
1009  if x and y are distinct eigenvectors. The matrix M is optional and is
1010  assumed to be the identity if left unset.
1011 
1012  The efficiency of LOBPCG relies on the availability of a suitable
1013  preconditioner for the matrix A. The preconditioner is supplied through the
1014  SetPreconditioner() method. It should be noted that the operator used with
1015  the preconditioner need not be A itself.
1016 
1017  For more information regarding LOBPCG see "Block Locally Optimal
1018  Preconditioned Eigenvalue Xolvers (BLOPEX) in Hypre and PETSc" by
1019  A. Knyazev, M. Argentati, I. Lashuk, and E. Ovtchinnikov, SISC, 29(5),
1020  2224-2239, 2007.
1021 */
1023 {
1024 private:
1025  MPI_Comm comm;
1026  int myid;
1027  int numProcs;
1028  int nev; // Number of desired eigenmodes
1029  int seed; // Random seed used for initial vectors
1030 
1031  HYPRE_Int glbSize; // Global number of DoFs in the linear system
1032  HYPRE_Int * part; // Row partitioning of the linear system
1033 
1034  // Pointer to HYPRE's solver struct
1035  HYPRE_Solver lobpcg_solver;
1036 
1037  // Interface for matrix storage type
1038  mv_InterfaceInterpreter interpreter;
1039 
1040  // Interface for setting up and performing matrix-vector products
1041  HYPRE_MatvecFunctions matvec_fn;
1042 
1043  // Eigenvalues
1044  Array<double> eigenvalues;
1045 
1046  // Forward declaration
1047  class HypreMultiVector;
1048 
1049  // MultiVector to store eigenvectors
1050  HypreMultiVector * multi_vec;
1051 
1052  // Empty vectors used to setup the matrices and preconditioner
1053  HypreParVector * x;
1054 
1055  // An optional operator which projects vectors into a desired subspace
1056  Operator * subSpaceProj;
1057 
1058  /// Internal class to represent a set of eigenvectors
1059  class HypreMultiVector
1060  {
1061  private:
1062  // Pointer to hypre's multi-vector object
1063  mv_MultiVectorPtr mv_ptr;
1064 
1065  // Wrappers for each member of the multivector
1066  HypreParVector ** hpv;
1067 
1068  // Number of vectors in the multivector
1069  int nv;
1070 
1071  public:
1072  HypreMultiVector(int n, HypreParVector & v,
1073  mv_InterfaceInterpreter & interpreter);
1074  ~HypreMultiVector();
1075 
1076  /// Set random values
1077  void Randomize(HYPRE_Int seed);
1078 
1079  /// Extract a single HypreParVector object
1080  HypreParVector & GetVector(unsigned int i);
1081 
1082  /// Transfers ownership of data to returned array of vectors
1083  HypreParVector ** StealVectors();
1084 
1085  operator mv_MultiVectorPtr() const { return mv_ptr; }
1086 
1087  mv_MultiVectorPtr & GetMultiVector() { return mv_ptr; }
1088  };
1089 
1090  static void * OperatorMatvecCreate( void *A, void *x );
1091  static HYPRE_Int OperatorMatvec( void *matvec_data,
1092  HYPRE_Complex alpha,
1093  void *A,
1094  void *x,
1095  HYPRE_Complex beta,
1096  void *y );
1097  static HYPRE_Int OperatorMatvecDestroy( void *matvec_data );
1098 
1099  static HYPRE_Int PrecondSolve(void *solver,
1100  void *A,
1101  void *b,
1102  void *x);
1103  static HYPRE_Int PrecondSetup(void *solver,
1104  void *A,
1105  void *b,
1106  void *x);
1107 
1108 public:
1109  HypreLOBPCG(MPI_Comm comm);
1110  ~HypreLOBPCG();
1111 
1112  void SetTol(double tol);
1113  void SetRelTol(double rel_tol);
1114  void SetMaxIter(int max_iter);
1115  void SetPrintLevel(int logging);
1116  void SetNumModes(int num_eigs) { nev = num_eigs; }
1117  void SetPrecondUsageMode(int pcg_mode);
1118  void SetRandomSeed(int s) { seed = s; }
1119  void SetInitialVectors(int num_vecs, HypreParVector ** vecs);
1120 
1121  // The following four methods support general operators
1122  void SetPreconditioner(Solver & precond);
1123  void SetOperator(Operator & A);
1124  void SetMassMatrix(Operator & M);
1125  void SetSubSpaceProjector(Operator & proj) { subSpaceProj = &proj; }
1126 
1127  /// Solve the eigenproblem
1128  void Solve();
1129 
1130  /// Collect the converged eigenvalues
1131  void GetEigenvalues(Array<double> & eigenvalues);
1132 
1133  /// Extract a single eigenvector
1134  HypreParVector & GetEigenvector(unsigned int i);
1135 
1136  /// Transfer ownership of the converged eigenvectors
1137  HypreParVector ** StealEigenvectors() { return multi_vec->StealVectors(); }
1138 };
1139 
1140 /** AME eigenvalue solver in hypre
1141 
1142  The Auxiliary space Maxwell Eigensolver (AME) is designed to find
1143  the lowest eigenmodes of the generalized eigenvalue problem:
1144  Curl Curl x = lambda M x
1145  where the Curl Curl operator is discretized using Nedelec finite element
1146  basis functions. Properties of this discretization are essential to
1147  eliminating the large null space of the Curl Curl operator.
1148 
1149  This eigensolver relies upon the LOBPCG eigensolver internally. It is also
1150  expected that the preconditioner supplied to this method will be the
1151  HypreAMS preconditioner defined above.
1152 
1153  As with LOBPCG, the operator set in the preconditioner need not be the same
1154  as A. This flexibility may be useful in solving eigenproblems which bare a
1155  strong resemblance to the Curl Curl problems for which AME is designed.
1156 
1157  Unlike LOBPCG, this eigensolver requires that the mass matrix be set.
1158  It is possible to circumvent this by passing an identity operator as the
1159  mass matrix but it seems unlikely that this would be useful so it is not the
1160  default behavior.
1161 */
1163 {
1164 private:
1165  int myid;
1166  int numProcs;
1167  int nev; // Number of desired eigenmodes
1168  bool setT;
1169 
1170  // Pointer to HYPRE's AME solver struct
1171  HYPRE_Solver ame_solver;
1172 
1173  // Pointer to HYPRE's AMS solver struct
1174  HypreSolver * ams_precond;
1175 
1176  // Eigenvalues
1177  HYPRE_Real * eigenvalues;
1178 
1179  // MultiVector to store eigenvectors
1180  HYPRE_ParVector * multi_vec;
1181 
1182  HypreParVector ** eigenvectors;
1183 
1184  void createDummyVectors();
1185 
1186 public:
1187  HypreAME(MPI_Comm comm);
1188  ~HypreAME();
1189 
1190  void SetTol(double tol);
1191  void SetRelTol(double rel_tol);
1192  void SetMaxIter(int max_iter);
1193  void SetPrintLevel(int logging);
1194  void SetNumModes(int num_eigs);
1195 
1196  // The following four methods support operators of type HypreParMatrix.
1197  void SetPreconditioner(HypreSolver & precond);
1198  void SetOperator(HypreParMatrix & A);
1199  void SetMassMatrix(HypreParMatrix & M);
1200 
1201  /// Solve the eigenproblem
1202  void Solve();
1203 
1204  /// Collect the converged eigenvalues
1205  void GetEigenvalues(Array<double> & eigenvalues);
1206 
1207  /// Extract a single eigenvector
1208  HypreParVector & GetEigenvector(unsigned int i);
1209 
1210  /// Transfer ownership of the converged eigenvectors
1212 };
1213 
1214 }
1215 
1216 #endif // MFEM_USE_MPI
1217 
1218 #endif
Type GetType() const
Definition: hypre.hpp:537
virtual ~HypreBoomerAMG()
Definition: hypre.cpp:2743
const HYPRE_Int * ColPart() const
Returns the column partitioning (const version)
Definition: hypre.hpp:380
void SetTol(double tol)
Definition: hypre.cpp:2182
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2348
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:1373
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:332
Vector * GlobalVector() const
Returns the global vector in each processor.
Definition: hypre.cpp:176
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:936
HypreParVector * X0
FIR Filter Temporary Vectors.
Definition: hypre.hpp:577
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:3412
double min_eig_est
Minimal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:607
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:800
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:968
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3193
void Print(const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0)
Prints the locally owned rows in parallel.
Definition: hypre.cpp:1396
MPI_Comm GetComm()
MPI communicator.
Definition: hypre.hpp:117
HYPRE_Int N() const
Returns the global number of columns.
Definition: hypre.hpp:384
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1563
virtual ~HypreEuclid()
Definition: hypre.cpp:2484
int setup_called
Was hypre&#39;s Setup function called already?
Definition: hypre.hpp:680
HYPRE_Int * ColPart()
Returns the column partitioning.
Definition: hypre.hpp:372
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to &#39;master&#39;.
Definition: hypre.cpp:822
void Read_IJMatrix(MPI_Comm comm, const char *fname)
Read a matrix saved as a HYPRE_IJMatrix.
Definition: hypre.cpp:1416
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:1086
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:1080
HypreParMatrix & Add(const double beta, const HypreParMatrix &B)
Definition: hypre.hpp:483
HypreParVector * B
Right-hand side and solution vector.
Definition: hypre.hpp:677
HypreDiagScale(HypreParMatrix &A)
Definition: hypre.hpp:811
double window_params[3]
Parameters for windowing function of FIR filter.
Definition: hypre.hpp:609
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: hypre.hpp:447
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:961
void SetInitialVectors(int num_vecs, HypreParVector **vecs)
Definition: hypre.cpp:3418
HypreADS(HypreParMatrix &A, ParFiniteElementSpace *face_fespace)
Definition: hypre.cpp:2953
void SetWindowByName(const char *window_name)
Convenience function for setting canonical windowing parameters.
Definition: hypre.cpp:1847
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:3679
void SetPreconditioner(HypreSolver &precond)
Definition: hypre.cpp:3645
virtual N_Vector ToNVector()
Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_PARHYP.
Definition: hypre.hpp:168
void SetTol(double tol)
Definition: hypre.cpp:3614
virtual ~HypreGMRES()
Definition: hypre.cpp:2424
void SetPreconditioner(Solver &precond)
Definition: hypre.cpp:3335
void SetMassMatrix(Operator &M)
Definition: hypre.cpp:3390
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetPrintLevel(int logging)
Definition: hypre.cpp:3320
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:283
virtual ~HypreAMS()
Definition: hypre.cpp:2933
void SetTol(double tol)
Definition: hypre.cpp:2323
void SetOperator(HypreParMatrix &A)
Definition: hypre.cpp:3651
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
GMRES Solve function.
Definition: hypre.hpp:782
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.hpp:1137
HypreParMatrix & operator+=(const HypreParMatrix &B)
Definition: hypre.hpp:476
bool pos_l1_norms
If set, take absolute values of the computed l1_norms.
Definition: hypre.hpp:603
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:1006
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2343
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
Definition: hypre.cpp:1183
double poly_fraction
Fraction of spectrum to smooth for polynomial relaxation.
Definition: hypre.hpp:591
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A * B.
Definition: hypre.cpp:1543
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2197
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:846
int poly_scale
Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}.
Definition: hypre.hpp:593
void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
Definition: hypre.hpp:462
void SetTol(double tol)
Definition: hypre.cpp:3298
HypreParVector * X1
Definition: hypre.hpp:577
The identity operator as a hypre solver.
Definition: hypre.hpp:793
void SetWindowParameters(double a, double b, double c)
Set parameters for windowing function for FIR smoother.
Definition: hypre.cpp:1862
void EliminateZeroRows()
If a row contains only zeros, set its diagonal to 1.
Definition: hypre.hpp:512
HYPRE_Int NNZ() const
Returns the global number of nonzeros.
Definition: hypre.hpp:364
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s GMRES.
Definition: hypre.cpp:2356
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
Definition: hypre.hpp:655
void Print(const char *fname) const
Prints the locally owned rows in parallel.
Definition: hypre.cpp:219
void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
Definition: hypre.cpp:1617
void SetSystemsOptions(int dim)
Definition: hypre.cpp:2625
HypreLOBPCG(MPI_Comm comm)
Definition: hypre.cpp:3268
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:3030
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
Definition: hypre.cpp:1986
hypre_ParVector * StealParVector()
Changes the ownership of the the vector.
Definition: hypre.hpp:134
The BoomerAMG solver in hypre.
Definition: hypre.hpp:873
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:417
void SetSymmetry(int sym)
Definition: hypre.cpp:2453
int dim
Definition: ex3.cpp:48
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:449
virtual ~HypreParaSails()
Definition: hypre.cpp:2458
The ParaSails preconditioner in hypre.
Definition: hypre.hpp:824
HYPRE_Int GetGlobalNumCols() const
Return the global number of columns.
Definition: hypre.hpp:421
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3314
double * l1_norms
l1 norms of the rows of A
Definition: hypre.hpp:601
Data type sparse matrix.
Definition: sparsemat.hpp:40
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:344
void SetLogging(int logging)
Definition: hypre.cpp:2192
void SetZeroInintialIterate()
non-hypre setting
Definition: hypre.hpp:728
Jacobi preconditioner in hypre.
Definition: hypre.hpp:807
virtual ~HypreSolver()
Definition: hypre.cpp:2164
HyprePCG(HypreParMatrix &_A)
Definition: hypre.cpp:2171
virtual HYPRE_PtrToParSolverFcn SolveFcn() const =0
hypre&#39;s internal Solve function
void SetRelTol(double rel_tol)
Definition: hypre.cpp:3620
void SetKDim(int dim)
Definition: hypre.cpp:2333
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:146
void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0)
Definition: hypre.cpp:2210
int to_int(const std::string &str)
Definition: text.hpp:70
void SetPrintLevel(int print_level)
Definition: hypre.hpp:914
void GetBlocks(Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
Definition: hypre.cpp:966
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
Definition: hypre.cpp:2110
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
PCG Solve function.
Definition: hypre.hpp:744
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.cpp:3726
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:988
HYPRE_Int Randomize(HYPRE_Int seed)
Set random values.
Definition: hypre.cpp:214
double relax_weight
Damping coefficient (usually &lt;= 1)
Definition: hypre.hpp:585
Parallel smoothers in hypre.
Definition: hypre.hpp:567
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:2707
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:924
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:920
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:798
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:571
HypreAMS(HypreParMatrix &A, ParFiniteElementSpace *edge_fespace)
Definition: hypre.cpp:2754
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:816
Dynamic 2D array using row-major layout.
Definition: array.hpp:304
virtual void SetOperator(const Operator &op)
Definition: hypre.cpp:1869
virtual HYPRE_PtrToParSolverFcn SetupFcn() const =0
hypre&#39;s internal Setup function
void SetLogging(int logging)
Definition: hypre.cpp:2338
virtual ~HypreADS()
Definition: hypre.cpp:3171
HypreParMatrix()
An empty matrix to be used as a reference to an existing matrix.
Definition: hypre.cpp:305
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:814
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2187
HypreParMatrix * GetData()
Definition: hypre.hpp:819
HYPRE_Int * RowPart()
Returns the row partitioning.
Definition: hypre.hpp:368
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
Definition: hypre.cpp:1827
void BooleanMult(int alpha, const int *x, int beta, int *y)
Definition: hypre.hpp:454
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
HypreParVector * X
Definition: hypre.hpp:677
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:992
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:135
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:3400
PCG solver in hypre.
Definition: hypre.hpp:706
GMRES solver in hypre.
Definition: hypre.hpp:755
signed char OwnsOffd() const
Get offd ownership flag.
Definition: hypre.hpp:350
const HYPRE_Int * RowPart() const
Returns the row partitioning (const version)
Definition: hypre.hpp:376
signed char OwnsColMap() const
Get colmap ownership flag.
Definition: hypre.hpp:352
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2948
HypreParVector * X
Definition: hypre.hpp:573
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:3673
void SetNumModes(int num_eigs)
Definition: hypre.cpp:3606
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:839
HypreParMatrix(hypre_ParCSRMatrix *a, bool owner=true)
Converts hypre&#39;s format to HypreParMatrix.
Definition: hypre.hpp:239
virtual ~HypreDiagScale()
Definition: hypre.hpp:820
int GetNumRows() const
Returns the number of rows in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:403
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.hpp:695
HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A + B.
Definition: hypre.cpp:1554
double max_eig_est
Maximal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:605
A class to initialize the size of a Tensor.
Definition: dtensor.hpp:55
void SetRelTol(double rel_tol)
Definition: hypre.cpp:3304
HypreAME(MPI_Comm comm)
Definition: hypre.cpp:3564
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:837
virtual ~HypreIdentity()
Definition: hypre.hpp:803
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:24
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3630
void SetFIRCoefficients(double max_eig)
Compute window and Chebyshev coefficients for given polynomial order.
Definition: hypre.cpp:1946
void SetRandomSeed(int s)
Definition: hypre.hpp:1118
signed char OwnsDiag() const
Get diag ownership flag.
Definition: hypre.hpp:348
double * fir_coeffs
Combined coefficients for windowing and Chebyshev polynomials.
Definition: hypre.hpp:612
void Threshold(double threshold=0.0)
Remove values smaller in absolute value than some threshold.
Definition: hypre.cpp:1297
double InnerProduct(HypreParVector *x, HypreParVector *y)
Definition: hypre.cpp:253
HYPRE_Int GlobalSize()
Returns the global number of rows.
Definition: hypre.hpp:125
virtual ~HypreSmoother()
Definition: hypre.cpp:2076
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:98
void SetSubSpaceProjector(Operator &proj)
Definition: hypre.hpp:1125
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2608
void SetData(double *_data)
Sets the data of the Vector and the hypre_ParVector to _data.
Definition: hypre.cpp:208
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
Definition: hypre.cpp:1218
void SetTaubinOptions(double lambda, double mu, int iter)
Set parameters for Taubin&#39;s lambda-mu method.
Definition: hypre.cpp:1839
HypreParVector * B
Right-hand side and solution vectors.
Definition: hypre.hpp:573
void SetMassMatrix(HypreParMatrix &M)
Definition: hypre.cpp:3666
HYPRE_Int M() const
Returns the global number of rows.
Definition: hypre.hpp:382
void SetZeroInintialIterate()
non-hypre setting
Definition: hypre.hpp:773
virtual ~HyprePCG()
Definition: hypre.cpp:2299
int GetNumCols() const
Returns the number of columns in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:410
void SetOperator(Operator &A)
Definition: hypre.cpp:3344
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2202
void SetPolyOptions(int poly_order, double poly_fraction)
Set parameters for polynomial smoothing.
Definition: hypre.cpp:1833
void GetNumIterations(int &num_iterations)
Definition: hypre.hpp:730
void CopyColStarts()
Definition: hypre.cpp:883
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:866
int relax_times
Number of relaxation sweeps.
Definition: hypre.hpp:583
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
Definition: hypre.cpp:1259
void SetPrintLevel(int logging)
Definition: hypre.cpp:3636
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:264
void SetOwnership(int own)
Sets ownership of the internal hypre_ParVector.
Definition: hypre.hpp:137
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:670
HypreParVector & operator=(double d)
Set constant values.
Definition: hypre.cpp:186
HypreEuclid(HypreParMatrix &A)
Definition: hypre.cpp:2464
const double alpha
Definition: ex15.cpp:339
Vector data type.
Definition: vector.hpp:48
void PrintCommPkg(std::ostream &out=mfem::out) const
Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
Definition: hypre.cpp:1437
HypreParMatrix & operator=(double value)
Initialize all entries with value.
Definition: hypre.hpp:469
HypreGMRES(HypreParMatrix &_A)
Definition: hypre.cpp:2305
hypre_ParCSRMatrix * StealData()
Changes the ownership of the the matrix.
Definition: hypre.cpp:832
void SetPrecondUsageMode(int pcg_mode)
Definition: hypre.cpp:3329
virtual ~HypreParMatrix()
Calls hypre&#39;s destroy function.
Definition: hypre.hpp:535
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:674
ID for class HypreParMatrix.
Definition: operator.hpp:139
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2328
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:960
HYPRE_Int * GetColStarts() const
Return the parallel column partitioning array.
Definition: hypre.hpp:432
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:3715
double omega
SOR parameter (usually in (0,2))
Definition: hypre.hpp:587
HYPRE_Int * GetRowStarts() const
Return the parallel row partitioning array.
Definition: hypre.hpp:427
Base class for solvers.
Definition: operator.hpp:279
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:959
int GetOwnership() const
Gets ownership of the internal hypre_ParVector.
Definition: hypre.hpp:140
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
PCG Setup function.
Definition: hypre.hpp:741
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
Abstract operator.
Definition: operator.hpp:21
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:864
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:922
~HypreParVector()
Calls hypre&#39;s destroy function.
Definition: hypre.cpp:224
void SetNumModes(int num_eigs)
Definition: hypre.hpp:1116
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2223
HypreParVector * Z
Definition: hypre.hpp:575
void SetSingularProblem()
Set this option when solving a curl-curl problem with zero mass term.
Definition: hypre.hpp:954
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:3457
void Read(MPI_Comm comm, const char *fname)
Reads the matrix from a file.
Definition: hypre.cpp:1401
HypreParVector * V
Temporary vectors.
Definition: hypre.hpp:575
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:1821
HypreParaSails(HypreParMatrix &A)
Definition: hypre.cpp:2430
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
GMRES Setup function.
Definition: hypre.hpp:779
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:25
int poly_order
Order of the smoothing polynomial.
Definition: hypre.hpp:589
HYPRE_Int * Partitioning()
Returns the parallel row/column partitioning.
Definition: hypre.hpp:122
double lambda
Taubin&#39;s lambda-mu method parameters.
Definition: hypre.hpp:596
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:994