MFEM  v4.1.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 #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  /** @brief Copy constructor for a ParCSR matrix which creates a deep copy of
329  structure and data from @a P. */
330  HypreParMatrix(const HypreParMatrix &P);
331 
332  /// Make this HypreParMatrix a reference to 'master'
333  void MakeRef(const HypreParMatrix &master);
334 
335  /// MPI communicator
336  MPI_Comm GetComm() const { return A->comm; }
337 
338  /// Typecasting to hypre's hypre_ParCSRMatrix*
339  operator hypre_ParCSRMatrix*() const { return A; }
340 #ifndef HYPRE_PAR_CSR_MATRIX_STRUCT
341  /// Typecasting to hypre's HYPRE_ParCSRMatrix, a.k.a. void *
342  operator HYPRE_ParCSRMatrix() { return (HYPRE_ParCSRMatrix) A; }
343 #endif
344  /// Changes the ownership of the the matrix
345  hypre_ParCSRMatrix* StealData();
346 
347  /// Explicitly set the three ownership flags, see docs for diagOwner etc.
348  void SetOwnerFlags(signed char diag, signed char offd, signed char colmap)
349  { diagOwner = diag, offdOwner = offd, colMapOwner = colmap; }
350 
351  /// Get diag ownership flag
352  signed char OwnsDiag() const { return diagOwner; }
353  /// Get offd ownership flag
354  signed char OwnsOffd() const { return offdOwner; }
355  /// Get colmap ownership flag
356  signed char OwnsColMap() const { return colMapOwner; }
357 
358  /** If the HypreParMatrix does not own the row-starts array, make a copy of
359  it that the HypreParMatrix will own. If the col-starts array is the same
360  as the row-starts array, col-starts is also replaced. */
361  void CopyRowStarts();
362  /** If the HypreParMatrix does not own the col-starts array, make a copy of
363  it that the HypreParMatrix will own. If the row-starts array is the same
364  as the col-starts array, row-starts is also replaced. */
365  void CopyColStarts();
366 
367  /// Returns the global number of nonzeros
368  inline HYPRE_Int NNZ() const { return A->num_nonzeros; }
369  /// Returns the row partitioning
370  /** See @ref hypre_partitioning_descr "here" for a description of the
371  partitioning array. */
372  inline HYPRE_Int *RowPart() { return A->row_starts; }
373  /// Returns the column partitioning
374  /** See @ref hypre_partitioning_descr "here" for a description of the
375  partitioning array. */
376  inline HYPRE_Int *ColPart() { return A->col_starts; }
377  /// Returns the row partitioning (const version)
378  /** See @ref hypre_partitioning_descr "here" for a description of the
379  partitioning array. */
380  inline const HYPRE_Int *RowPart() const { return A->row_starts; }
381  /// Returns the column partitioning (const version)
382  /** See @ref hypre_partitioning_descr "here" for a description of the
383  partitioning array. */
384  inline const HYPRE_Int *ColPart() const { return A->col_starts; }
385  /// Returns the global number of rows
386  inline HYPRE_Int M() const { return A->global_num_rows; }
387  /// Returns the global number of columns
388  inline HYPRE_Int N() const { return A->global_num_cols; }
389 
390  /// Get the local diagonal of the matrix.
391  void GetDiag(Vector &diag) const;
392  /// Get the local diagonal block. NOTE: 'diag' will not own any data.
393  void GetDiag(SparseMatrix &diag) const;
394  /// Get the local off-diagonal block. NOTE: 'offd' will not own any data.
395  void GetOffd(SparseMatrix &offd, HYPRE_Int* &cmap) const;
396 
397  /** Split the matrix into M x N equally sized blocks of parallel matrices.
398  The size of 'blocks' must already be set to M x N. */
399  void GetBlocks(Array2D<HypreParMatrix*> &blocks,
400  bool interleaved_rows = false,
401  bool interleaved_cols = false) const;
402 
403  /// Returns the transpose of *this
404  HypreParMatrix * Transpose() const;
405 
406  /// Returns the number of rows in the diagonal block of the ParCSRMatrix
407  int GetNumRows() const
408  {
409  return internal::to_int(
410  hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(A)));
411  }
412 
413  /// Returns the number of columns in the diagonal block of the ParCSRMatrix
414  int GetNumCols() const
415  {
416  return internal::to_int(
417  hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(A)));
418  }
419 
420  /// Return the global number of rows
421  HYPRE_Int GetGlobalNumRows() const
422  { return hypre_ParCSRMatrixGlobalNumRows(A); }
423 
424  /// Return the global number of columns
425  HYPRE_Int GetGlobalNumCols() const
426  { return hypre_ParCSRMatrixGlobalNumCols(A); }
427 
428  /// Return the parallel row partitioning array.
429  /** See @ref hypre_partitioning_descr "here" for a description of the
430  partitioning array. */
431  HYPRE_Int *GetRowStarts() const { return hypre_ParCSRMatrixRowStarts(A); }
432 
433  /// Return the parallel column partitioning array.
434  /** See @ref hypre_partitioning_descr "here" for a description of the
435  partitioning array. */
436  HYPRE_Int *GetColStarts() const { return hypre_ParCSRMatrixColStarts(A); }
437 
438  /// Computes y = alpha * A * x + beta * y
439  HYPRE_Int Mult(HypreParVector &x, HypreParVector &y,
440  double alpha = 1.0, double beta = 0.0);
441  /// Computes y = alpha * A * x + beta * y
442  HYPRE_Int Mult(HYPRE_ParVector x, HYPRE_ParVector y,
443  double alpha = 1.0, double beta = 0.0);
444  /// Computes y = alpha * A^t * x + beta * y
446  double alpha = 1.0, double beta = 0.0);
447 
448  void Mult(double a, const Vector &x, double b, Vector &y) const;
449  void MultTranspose(double a, const Vector &x, double b, Vector &y) const;
450 
451  virtual void Mult(const Vector &x, Vector &y) const
452  { Mult(1.0, x, 0.0, y); }
453  virtual void MultTranspose(const Vector &x, Vector &y) const
454  { MultTranspose(1.0, x, 0.0, y); }
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 /** Eliminate essential BC specified by 'ess_dof_list' from the solution X to
574  the r.h.s. B. Here A is a matrix with eliminated BC, while Ae is such that
575  (A+Ae) is the original (Neumann) matrix before elimination. */
576 void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae,
577  const Array<int> &ess_dof_list, const Vector &X, Vector &B);
578 
579 
580 /// Parallel smoothers in hypre
581 class HypreSmoother : public Solver
582 {
583 protected:
584  /// The linear system matrix
586  /// Right-hand side and solution vectors
587  mutable HypreParVector *B, *X;
588  /// Temporary vectors
589  mutable HypreParVector *V, *Z;
590  /// FIR Filter Temporary Vectors
591  mutable HypreParVector *X0, *X1;
592 
593  /** Smoother type from hypre_ParCSRRelax() in ams.c plus extensions, see the
594  enumeration Type below. */
595  int type;
596  /// Number of relaxation sweeps
598  /// Damping coefficient (usually <= 1)
599  double relax_weight;
600  /// SOR parameter (usually in (0,2))
601  double omega;
602  /// Order of the smoothing polynomial
604  /// Fraction of spectrum to smooth for polynomial relaxation
606  /// Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}
608 
609  /// Taubin's lambda-mu method parameters
610  double lambda;
611  double mu;
613 
614  /// l1 norms of the rows of A
615  double *l1_norms;
616  /// If set, take absolute values of the computed l1_norms
618  /// Maximal eigenvalue estimate for polynomial smoothing
619  double max_eig_est;
620  /// Minimal eigenvalue estimate for polynomial smoothing
621  double min_eig_est;
622  /// Parameters for windowing function of FIR filter
623  double window_params[3];
624 
625  /// Combined coefficients for windowing and Chebyshev polynomials.
626  double* fir_coeffs;
627 
628 public:
629  /** Hypre smoother types:
630  0 = Jacobi
631  1 = l1-scaled Jacobi
632  2 = l1-scaled block Gauss-Seidel/SSOR
633  4 = truncated l1-scaled block Gauss-Seidel/SSOR
634  5 = lumped Jacobi
635  6 = Gauss-Seidel
636  16 = Chebyshev
637  1001 = Taubin polynomial smoother
638  1002 = FIR polynomial smoother. */
639  enum Type { Jacobi = 0, l1Jacobi = 1, l1GS = 2, l1GStr = 4, lumpedJacobi = 5,
640  GS = 6, Chebyshev = 16, Taubin = 1001, FIR = 1002
641  };
642 
643  HypreSmoother();
644 
646  int relax_times = 1, double relax_weight = 1.0,
647  double omega = 1.0, int poly_order = 2,
648  double poly_fraction = .3);
649 
650  /// Set the relaxation type and number of sweeps
652  /// Set SOR-related parameters
653  void SetSOROptions(double relax_weight, double omega);
654  /// Set parameters for polynomial smoothing
655  void SetPolyOptions(int poly_order, double poly_fraction);
656  /// Set parameters for Taubin's lambda-mu method
657  void SetTaubinOptions(double lambda, double mu, int iter);
658 
659  /// Convenience function for setting canonical windowing parameters
660  void SetWindowByName(const char* window_name);
661  /// Set parameters for windowing function for FIR smoother.
662  void SetWindowParameters(double a, double b, double c);
663  /// Compute window and Chebyshev coefficients for given polynomial order.
664  void SetFIRCoefficients(double max_eig);
665 
666  /// After computing l1-norms, replace them with their absolute values.
667  /** By default, the l1-norms take their sign from the corresponding diagonal
668  entries in the associated matrix. */
669  void SetPositiveDiagonal(bool pos = true) { pos_l1_norms = pos; }
670 
671  /** Set/update the associated operator. Must be called after setting the
672  HypreSmoother type and options. */
673  virtual void SetOperator(const Operator &op);
674 
675  /// Relax the linear system Ax=b
676  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
677  virtual void Mult(const Vector &b, Vector &x) const;
678 
679  virtual ~HypreSmoother();
680 };
681 
682 
683 /// Abstract class for hypre's solvers and preconditioners
684 class HypreSolver : public Solver
685 {
686 public:
687  /// How to treat errors returned by hypre function calls.
689  {
690  IGNORE_HYPRE_ERRORS, ///< Ignore hypre errors (see e.g. HypreADS)
691  WARN_HYPRE_ERRORS, ///< Issue warnings on hypre errors
692  ABORT_HYPRE_ERRORS ///< Abort on hypre errors (default in base class)
693  };
694 
695 protected:
696  /// The linear system matrix
698 
699  /// Right-hand side and solution vector
700  mutable HypreParVector *B, *X;
701 
702  /// Was hypre's Setup function called already?
703  mutable int setup_called;
704 
705  /// How to treat hypre errors.
707 
708 public:
709  HypreSolver();
710 
712 
713  /// Typecast to HYPRE_Solver -- return the solver
714  virtual operator HYPRE_Solver() const = 0;
715 
716  /// hypre's internal Setup function
717  virtual HYPRE_PtrToParSolverFcn SetupFcn() const = 0;
718  /// hypre's internal Solve function
719  virtual HYPRE_PtrToParSolverFcn SolveFcn() const = 0;
720 
721  virtual void SetOperator(const Operator &op)
722  { mfem_error("HypreSolvers do not support SetOperator!"); }
723 
724  /// Solve the linear system Ax=b
725  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
726  virtual void Mult(const Vector &b, Vector &x) const;
727 
728  /** @brief Set the behavior for treating hypre errors, see the ErrorMode
729  enum. The default mode in the base class is ABORT_HYPRE_ERRORS. */
730  /** Currently, there are three cases in derived classes where the error flag
731  is set to IGNORE_HYPRE_ERRORS:
732  * in the method HypreBoomerAMG::SetElasticityOptions(), and
733  * in the constructor of classes HypreAMS and HypreADS.
734  The reason for this is that a nonzero hypre error is returned) when
735  hypre_ParCSRComputeL1Norms() encounters zero row in a matrix, which is
736  expected in some cases with the above solvers. */
737  void SetErrorMode(ErrorMode err_mode) const { error_mode = err_mode; }
738 
739  virtual ~HypreSolver();
740 };
741 
742 /// PCG solver in hypre
743 class HyprePCG : public HypreSolver
744 {
745 private:
746  HYPRE_Solver pcg_solver;
747 
748  HypreSolver * precond;
749 
750 public:
751  HyprePCG(MPI_Comm comm);
752 
754 
755  virtual void SetOperator(const Operator &op);
756 
757  void SetTol(double tol);
758  void SetMaxIter(int max_iter);
759  void SetLogging(int logging);
760  void SetPrintLevel(int print_lvl);
761 
762  /// Set the hypre solver to be used as a preconditioner
763  void SetPreconditioner(HypreSolver &precond);
764 
765  /** Use the L2 norm of the residual for measuring PCG convergence, plus
766  (optionally) 1) periodically recompute true residuals from scratch; and
767  2) enable residual-based stopping criteria. */
768  void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0);
769 
770  /// non-hypre setting
772 
773  void GetNumIterations(int &num_iterations)
774  {
775  HYPRE_Int num_it;
776  HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_it);
777  num_iterations = internal::to_int(num_it);
778  }
779 
780  /// The typecast to HYPRE_Solver returns the internal pcg_solver
781  virtual operator HYPRE_Solver() const { return pcg_solver; }
782 
783  /// PCG Setup function
784  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
785  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSetup; }
786  /// PCG Solve function
787  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
788  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSolve; }
789 
790  /// Solve Ax=b with hypre's PCG
791  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
792  using HypreSolver::Mult;
793 
794  virtual ~HyprePCG();
795 };
796 
797 /// GMRES solver in hypre
798 class HypreGMRES : public HypreSolver
799 {
800 private:
801  HYPRE_Solver gmres_solver;
802 
803  HypreSolver * precond;
804 
805  /// Default, generally robust, GMRES options
806  void SetDefaultOptions();
807 
808 public:
809  HypreGMRES(MPI_Comm comm);
810 
812 
813  virtual void SetOperator(const Operator &op);
814 
815  void SetTol(double tol);
816  void SetMaxIter(int max_iter);
817  void SetKDim(int dim);
818  void SetLogging(int logging);
819  void SetPrintLevel(int print_lvl);
820 
821  /// Set the hypre solver to be used as a preconditioner
822  void SetPreconditioner(HypreSolver &precond);
823 
824  /// non-hypre setting
826 
827  /// The typecast to HYPRE_Solver returns the internal gmres_solver
828  virtual operator HYPRE_Solver() const { return gmres_solver; }
829 
830  /// GMRES Setup function
831  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
832  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSetup; }
833  /// GMRES Solve function
834  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
835  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSolve; }
836 
837  /// Solve Ax=b with hypre's GMRES
838  virtual void Mult (const HypreParVector &b, HypreParVector &x) const;
839  using HypreSolver::Mult;
840 
841  virtual ~HypreGMRES();
842 };
843 
844 /// The identity operator as a hypre solver
846 {
847 public:
848  virtual operator HYPRE_Solver() const { return NULL; }
849 
850  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
851  { return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentitySetup; }
852  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
853  { return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentity; }
854 
855  virtual ~HypreIdentity() { }
856 };
857 
858 /// Jacobi preconditioner in hypre
860 {
861 public:
864  virtual operator HYPRE_Solver() const { return NULL; }
865 
866  virtual void SetOperator(const Operator &op);
867 
868  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
869  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScaleSetup; }
870  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
871  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScale; }
872 
873  HypreParMatrix* GetData() { return A; }
874  virtual ~HypreDiagScale() { }
875 };
876 
877 /// The ParaSails preconditioner in hypre
879 {
880 private:
881  HYPRE_Solver sai_precond;
882 
883  /// Default, generally robust, ParaSails options
884  void SetDefaultOptions();
885 
886  // If sai_precond is NULL, this method allocates it and sets default options.
887  // Otherwise the method saves the options from sai_precond, destroys it,
888  // allocates a new object, and sets its options to the saved values.
889  void ResetSAIPrecond(MPI_Comm comm);
890 
891 public:
892  HypreParaSails(MPI_Comm comm);
893 
895 
896  virtual void SetOperator(const Operator &op);
897 
898  void SetSymmetry(int sym);
899 
900  /// The typecast to HYPRE_Solver returns the internal sai_precond
901  virtual operator HYPRE_Solver() const { return sai_precond; }
902 
903  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
904  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSetup; }
905  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
906  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSolve; }
907 
908  virtual ~HypreParaSails();
909 };
910 
911 /** The Euclid preconditioner in Hypre
912 
913  Euclid implements the Parallel Incomplete LU factorization technique. For
914  more information see:
915 
916  "A Scalable Parallel Algorithm for Incomplete Factor Preconditioning" by
917  David Hysom and Alex Pothen, https://doi.org/10.1137/S1064827500376193
918 */
919 class HypreEuclid : public HypreSolver
920 {
921 private:
922  HYPRE_Solver euc_precond;
923 
924  /// Default, generally robust, Euclid options
925  void SetDefaultOptions();
926 
927  // If euc_precond is NULL, this method allocates it and sets default options.
928  // Otherwise the method saves the options from euc_precond, destroys it,
929  // allocates a new object, and sets its options to the saved values.
930  void ResetEuclidPrecond(MPI_Comm comm);
931 
932 public:
933  HypreEuclid(MPI_Comm comm);
934 
936 
937  virtual void SetOperator(const Operator &op);
938 
939  /// The typecast to HYPRE_Solver returns the internal euc_precond
940  virtual operator HYPRE_Solver() const { return euc_precond; }
941 
942  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
943  { return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSetup; }
944  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
945  { return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSolve; }
946 
947  virtual ~HypreEuclid();
948 };
949 
950 /// The BoomerAMG solver in hypre
952 {
953 private:
954  HYPRE_Solver amg_precond;
955 
956  /// Rigid body modes
958 
959  /// Finite element space for elasticity problems, see SetElasticityOptions()
960  ParFiniteElementSpace *fespace;
961 
962  /// Recompute the rigid-body modes vectors (in the rbms array)
963  void RecomputeRBMs();
964 
965  /// Default, generally robust, BoomerAMG options
966  void SetDefaultOptions();
967 
968  // If amg_precond is NULL, allocates it and sets default options.
969  // Otherwise saves the options from amg_precond, destroys it, allocates a new
970  // one, and sets its options to the saved values.
971  void ResetAMGPrecond();
972 
973 public:
974  HypreBoomerAMG();
975 
977 
978  virtual void SetOperator(const Operator &op);
979 
980  /** More robust options for systems, such as elasticity. Note that BoomerAMG
981  assumes Ordering::byVDIM in the finite element space used to generate the
982  matrix A. */
983  void SetSystemsOptions(int dim);
984 
985  /** A special elasticity version of BoomerAMG that takes advantage of
986  geometric rigid body modes and could perform better on some problems, see
987  "Improving algebraic multigrid interpolation operators for linear
988  elasticity problems", Baker, Kolev, Yang, NLAA 2009, DOI:10.1002/nla.688.
989  As with SetSystemsOptions(), this solver assumes Ordering::byVDIM. */
991 
992  void SetPrintLevel(int print_level)
993  { HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level); }
994 
995  /// The typecast to HYPRE_Solver returns the internal amg_precond
996  virtual operator HYPRE_Solver() const { return amg_precond; }
997 
998  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
999  { return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSetup; }
1000  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1001  { return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSolve; }
1002 
1003  virtual ~HypreBoomerAMG();
1004 };
1005 
1006 /// Compute the discrete gradient matrix between the nodal linear and ND1 spaces
1007 HypreParMatrix* DiscreteGrad(ParFiniteElementSpace *edge_fespace,
1008  ParFiniteElementSpace *vert_fespace);
1009 /// Compute the discrete curl matrix between the ND1 and RT0 spaces
1010 HypreParMatrix* DiscreteCurl(ParFiniteElementSpace *face_fespace,
1011  ParFiniteElementSpace *edge_fespace);
1012 
1013 /// The Auxiliary-space Maxwell Solver in hypre
1014 class HypreAMS : public HypreSolver
1015 {
1016 private:
1017  /// Constuct AMS solver from finite element space
1018  void Init(ParFiniteElementSpace *edge_space);
1019 
1020  HYPRE_Solver ams;
1021 
1022  /// Vertex coordinates
1023  HypreParVector *x, *y, *z;
1024  /// Discrete gradient matrix
1025  HypreParMatrix *G;
1026  /// Nedelec interpolation matrix and its components
1027  HypreParMatrix *Pi, *Pix, *Piy, *Piz;
1028 
1029 public:
1030  HypreAMS(ParFiniteElementSpace *edge_fespace);
1031 
1033 
1034  virtual void SetOperator(const Operator &op);
1035 
1036  void SetPrintLevel(int print_lvl);
1037 
1038  /// Set this option when solving a curl-curl problem with zero mass term
1039  void SetSingularProblem() { HYPRE_AMSSetBetaPoissonMatrix(ams, NULL); }
1040 
1041  /// The typecast to HYPRE_Solver returns the internal ams object
1042  virtual operator HYPRE_Solver() const { return ams; }
1043 
1044  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1045  { return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSetup; }
1046  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1047  { return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSolve; }
1048 
1049  virtual ~HypreAMS();
1050 };
1051 
1052 /// The Auxiliary-space Divergence Solver in hypre
1053 class HypreADS : public HypreSolver
1054 {
1055 private:
1056  /// Constuct ADS solver from finite element space
1057  void Init(ParFiniteElementSpace *face_fespace);
1058 
1059  HYPRE_Solver ads;
1060 
1061  /// Vertex coordinates
1062  HypreParVector *x, *y, *z;
1063  /// Discrete gradient matrix
1064  HypreParMatrix *G;
1065  /// Discrete curl matrix
1066  HypreParMatrix *C;
1067  /// Nedelec interpolation matrix and its components
1068  HypreParMatrix *ND_Pi, *ND_Pix, *ND_Piy, *ND_Piz;
1069  /// Raviart-Thomas interpolation matrix and its components
1070  HypreParMatrix *RT_Pi, *RT_Pix, *RT_Piy, *RT_Piz;
1071 
1072 public:
1073  HypreADS(ParFiniteElementSpace *face_fespace);
1074 
1076 
1077  virtual void SetOperator(const Operator &op);
1078 
1079  void SetPrintLevel(int print_lvl);
1080 
1081  /// The typecast to HYPRE_Solver returns the internal ads object
1082  virtual operator HYPRE_Solver() const { return ads; }
1083 
1084  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1085  { return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSetup; }
1086  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1087  { return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSolve; }
1088 
1089  virtual ~HypreADS();
1090 };
1091 
1092 /** LOBPCG eigenvalue solver in hypre
1093 
1094  The Locally Optimal Block Preconditioned Conjugate Gradient (LOBPCG)
1095  eigenvalue solver is designed to find the lowest eigenmodes of the
1096  generalized eigenvalue problem:
1097  A x = lambda M x
1098  where A is symmetric, potentially indefinite and M is symmetric positive
1099  definite. The eigenvectors are M-orthonormal, meaning that
1100  x^T M x = 1 and x^T M y = 0,
1101  if x and y are distinct eigenvectors. The matrix M is optional and is
1102  assumed to be the identity if left unset.
1103 
1104  The efficiency of LOBPCG relies on the availability of a suitable
1105  preconditioner for the matrix A. The preconditioner is supplied through the
1106  SetPreconditioner() method. It should be noted that the operator used with
1107  the preconditioner need not be A itself.
1108 
1109  For more information regarding LOBPCG see "Block Locally Optimal
1110  Preconditioned Eigenvalue Xolvers (BLOPEX) in Hypre and PETSc" by
1111  A. Knyazev, M. Argentati, I. Lashuk, and E. Ovtchinnikov, SISC, 29(5),
1112  2224-2239, 2007.
1113 */
1115 {
1116 private:
1117  MPI_Comm comm;
1118  int myid;
1119  int numProcs;
1120  int nev; // Number of desired eigenmodes
1121  int seed; // Random seed used for initial vectors
1122 
1123  HYPRE_Int glbSize; // Global number of DoFs in the linear system
1124  HYPRE_Int * part; // Row partitioning of the linear system
1125 
1126  // Pointer to HYPRE's solver struct
1127  HYPRE_Solver lobpcg_solver;
1128 
1129  // Interface for matrix storage type
1130  mv_InterfaceInterpreter interpreter;
1131 
1132  // Interface for setting up and performing matrix-vector products
1133  HYPRE_MatvecFunctions matvec_fn;
1134 
1135  // Eigenvalues
1136  Array<double> eigenvalues;
1137 
1138  // Forward declaration
1139  class HypreMultiVector;
1140 
1141  // MultiVector to store eigenvectors
1142  HypreMultiVector * multi_vec;
1143 
1144  // Empty vectors used to setup the matrices and preconditioner
1145  HypreParVector * x;
1146 
1147  // An optional operator which projects vectors into a desired subspace
1148  Operator * subSpaceProj;
1149 
1150  /// Internal class to represent a set of eigenvectors
1151  class HypreMultiVector
1152  {
1153  private:
1154  // Pointer to hypre's multi-vector object
1155  mv_MultiVectorPtr mv_ptr;
1156 
1157  // Wrappers for each member of the multivector
1158  HypreParVector ** hpv;
1159 
1160  // Number of vectors in the multivector
1161  int nv;
1162 
1163  public:
1164  HypreMultiVector(int n, HypreParVector & v,
1165  mv_InterfaceInterpreter & interpreter);
1166  ~HypreMultiVector();
1167 
1168  /// Set random values
1169  void Randomize(HYPRE_Int seed);
1170 
1171  /// Extract a single HypreParVector object
1172  HypreParVector & GetVector(unsigned int i);
1173 
1174  /// Transfers ownership of data to returned array of vectors
1175  HypreParVector ** StealVectors();
1176 
1177  operator mv_MultiVectorPtr() const { return mv_ptr; }
1178 
1179  mv_MultiVectorPtr & GetMultiVector() { return mv_ptr; }
1180  };
1181 
1182  static void * OperatorMatvecCreate( void *A, void *x );
1183  static HYPRE_Int OperatorMatvec( void *matvec_data,
1184  HYPRE_Complex alpha,
1185  void *A,
1186  void *x,
1187  HYPRE_Complex beta,
1188  void *y );
1189  static HYPRE_Int OperatorMatvecDestroy( void *matvec_data );
1190 
1191  static HYPRE_Int PrecondSolve(void *solver,
1192  void *A,
1193  void *b,
1194  void *x);
1195  static HYPRE_Int PrecondSetup(void *solver,
1196  void *A,
1197  void *b,
1198  void *x);
1199 
1200 public:
1201  HypreLOBPCG(MPI_Comm comm);
1202  ~HypreLOBPCG();
1203 
1204  void SetTol(double tol);
1205  void SetRelTol(double rel_tol);
1206  void SetMaxIter(int max_iter);
1207  void SetPrintLevel(int logging);
1208  void SetNumModes(int num_eigs) { nev = num_eigs; }
1209  void SetPrecondUsageMode(int pcg_mode);
1210  void SetRandomSeed(int s) { seed = s; }
1211  void SetInitialVectors(int num_vecs, HypreParVector ** vecs);
1212 
1213  // The following four methods support general operators
1214  void SetPreconditioner(Solver & precond);
1215  void SetOperator(Operator & A);
1216  void SetMassMatrix(Operator & M);
1217  void SetSubSpaceProjector(Operator & proj) { subSpaceProj = &proj; }
1218 
1219  /// Solve the eigenproblem
1220  void Solve();
1221 
1222  /// Collect the converged eigenvalues
1223  void GetEigenvalues(Array<double> & eigenvalues);
1224 
1225  /// Extract a single eigenvector
1226  HypreParVector & GetEigenvector(unsigned int i);
1227 
1228  /// Transfer ownership of the converged eigenvectors
1229  HypreParVector ** StealEigenvectors() { return multi_vec->StealVectors(); }
1230 };
1231 
1232 /** AME eigenvalue solver in hypre
1233 
1234  The Auxiliary space Maxwell Eigensolver (AME) is designed to find
1235  the lowest eigenmodes of the generalized eigenvalue problem:
1236  Curl Curl x = lambda M x
1237  where the Curl Curl operator is discretized using Nedelec finite element
1238  basis functions. Properties of this discretization are essential to
1239  eliminating the large null space of the Curl Curl operator.
1240 
1241  This eigensolver relies upon the LOBPCG eigensolver internally. It is also
1242  expected that the preconditioner supplied to this method will be the
1243  HypreAMS preconditioner defined above.
1244 
1245  As with LOBPCG, the operator set in the preconditioner need not be the same
1246  as A. This flexibility may be useful in solving eigenproblems which bare a
1247  strong resemblance to the Curl Curl problems for which AME is designed.
1248 
1249  Unlike LOBPCG, this eigensolver requires that the mass matrix be set.
1250  It is possible to circumvent this by passing an identity operator as the
1251  mass matrix but it seems unlikely that this would be useful so it is not the
1252  default behavior.
1253 */
1255 {
1256 private:
1257  int myid;
1258  int numProcs;
1259  int nev; // Number of desired eigenmodes
1260  bool setT;
1261 
1262  // Pointer to HYPRE's AME solver struct
1263  HYPRE_Solver ame_solver;
1264 
1265  // Pointer to HYPRE's AMS solver struct
1266  HypreSolver * ams_precond;
1267 
1268  // Eigenvalues
1269  HYPRE_Real * eigenvalues;
1270 
1271  // MultiVector to store eigenvectors
1272  HYPRE_ParVector * multi_vec;
1273 
1274  HypreParVector ** eigenvectors;
1275 
1276  void createDummyVectors();
1277 
1278 public:
1279  HypreAME(MPI_Comm comm);
1280  ~HypreAME();
1281 
1282  void SetTol(double tol);
1283  void SetRelTol(double rel_tol);
1284  void SetMaxIter(int max_iter);
1285  void SetPrintLevel(int logging);
1286  void SetNumModes(int num_eigs);
1287 
1288  // The following four methods support operators of type HypreParMatrix.
1289  void SetPreconditioner(HypreSolver & precond);
1290  void SetOperator(HypreParMatrix & A);
1291  void SetMassMatrix(HypreParMatrix & M);
1292 
1293  /// Solve the eigenproblem
1294  void Solve();
1295 
1296  /// Collect the converged eigenvalues
1297  void GetEigenvalues(Array<double> & eigenvalues);
1298 
1299  /// Extract a single eigenvector
1300  HypreParVector & GetEigenvector(unsigned int i);
1301 
1302  /// Transfer ownership of the converged eigenvectors
1304 };
1305 
1306 }
1307 
1308 #endif // MFEM_USE_MPI
1309 
1310 #endif
Type GetType() const
Definition: hypre.hpp:549
virtual ~HypreBoomerAMG()
Definition: hypre.cpp:3031
const HYPRE_Int * ColPart() const
Returns the column partitioning (const version)
Definition: hypre.hpp:384
void SetTol(double tol)
Definition: hypre.cpp:2305
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2506
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:1390
HypreADS(ParFiniteElementSpace *face_fespace)
Definition: hypre.cpp:3272
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:336
HypreEuclid(MPI_Comm comm)
Definition: hypre.cpp:2703
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:737
Vector * GlobalVector() const
Returns the global vector in each processor.
Definition: hypre.cpp:176
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1014
HypreParVector * X0
FIR Filter Temporary Vectors.
Definition: hypre.hpp:591
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:3763
double min_eig_est
Minimal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:621
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:852
ErrorMode
How to treat errors returned by hypre function calls.
Definition: hypre.hpp:688
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:1053
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3544
void Print(const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0)
Prints the locally owned rows in parallel.
Definition: hypre.cpp:1436
MPI_Comm GetComm()
MPI communicator.
Definition: hypre.hpp:117
HYPRE_Int N() const
Returns the global number of columns.
Definition: hypre.hpp:388
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1636
virtual ~HypreEuclid()
Definition: hypre.cpp:2766
int setup_called
Was hypre&#39;s Setup function called already?
Definition: hypre.hpp:703
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2285
HYPRE_Int * ColPart()
Returns the column partitioning.
Definition: hypre.hpp:376
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to &#39;master&#39;.
Definition: hypre.cpp:837
void Read_IJMatrix(MPI_Comm comm, const char *fname)
Read a matrix saved as a HYPRE_IJMatrix.
Definition: hypre.cpp:1456
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:1103
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:1097
HypreParMatrix & Add(const double beta, const HypreParMatrix &B)
Definition: hypre.hpp:487
HypreParVector * B
Right-hand side and solution vector.
Definition: hypre.hpp:700
HypreDiagScale(HypreParMatrix &A)
Definition: hypre.hpp:863
double window_params[3]
Parameters for windowing function of FIR filter.
Definition: hypre.hpp:623
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: hypre.hpp:451
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1046
void SetInitialVectors(int num_vecs, HypreParVector **vecs)
Definition: hypre.cpp:3769
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2461
void SetWindowByName(const char *window_name)
Convenience function for setting canonical windowing parameters.
Definition: hypre.cpp:1920
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:4030
Issue warnings on hypre errors.
Definition: hypre.hpp:691
void SetPreconditioner(HypreSolver &precond)
Definition: hypre.cpp:3996
HyprePCG(MPI_Comm comm)
Definition: hypre.cpp:2267
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:3965
HypreAMS(ParFiniteElementSpace *edge_fespace)
Definition: hypre.cpp:3041
virtual ~HypreGMRES()
Definition: hypre.cpp:2584
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:3506
void SetPreconditioner(Solver &precond)
Definition: hypre.cpp:3686
void SetMassMatrix(Operator &M)
Definition: hypre.cpp:3741
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetPrintLevel(int logging)
Definition: hypre.cpp:3671
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:504
virtual ~HypreAMS()
Definition: hypre.cpp:3252
void SetTol(double tol)
Definition: hypre.cpp:2481
void SetOperator(HypreParMatrix &A)
Definition: hypre.cpp:4002
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
GMRES Solve function.
Definition: hypre.hpp:834
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.hpp:1229
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:617
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
Definition: hypre.cpp:1610
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:1021
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2501
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
Definition: hypre.cpp:1200
double poly_fraction
Fraction of spectrum to smooth for polynomial relaxation.
Definition: hypre.hpp:605
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2320
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:861
int poly_scale
Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}.
Definition: hypre.hpp:607
void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
Definition: hypre.hpp:466
void SetTol(double tol)
Definition: hypre.cpp:3649
HypreParVector * X1
Definition: hypre.hpp:591
The identity operator as a hypre solver.
Definition: hypre.hpp:845
void SetWindowParameters(double a, double b, double c)
Set parameters for windowing function for FIR smoother.
Definition: hypre.cpp:1935
void EliminateZeroRows()
If a row contains only zeros, set its diagonal to 1.
Definition: hypre.hpp:516
HypreGMRES(MPI_Comm comm)
Definition: hypre.cpp:2430
HYPRE_Int NNZ() const
Returns the global number of nonzeros.
Definition: hypre.hpp:368
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2670
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s GMRES.
Definition: hypre.cpp:2516
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
Definition: hypre.hpp:669
void Print(const char *fname) const
Prints the locally owned rows in parallel.
Definition: hypre.cpp:216
void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
Definition: hypre.cpp:1690
void SetSystemsOptions(int dim)
Definition: hypre.cpp:2907
HypreLOBPCG(MPI_Comm comm)
Definition: hypre.cpp:3619
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:1926
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
Definition: hypre.cpp:2059
hypre_ParVector * StealParVector()
Changes the ownership of the the vector.
Definition: hypre.hpp:134
The BoomerAMG solver in hypre.
Definition: hypre.hpp:951
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2590
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:421
void SetSymmetry(int sym)
Definition: hypre.cpp:2692
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:453
virtual ~HypreParaSails()
Definition: hypre.cpp:2697
The ParaSails preconditioner in hypre.
Definition: hypre.hpp:878
HYPRE_Int GetGlobalNumCols() const
Return the global number of columns.
Definition: hypre.hpp:425
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3665
double * l1_norms
l1 norms of the rows of A
Definition: hypre.hpp:615
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:348
void SetLogging(int logging)
Definition: hypre.cpp:2315
void SetZeroInintialIterate()
non-hypre setting
Definition: hypre.hpp:771
Jacobi preconditioner in hypre.
Definition: hypre.hpp:859
virtual ~HypreSolver()
Definition: hypre.cpp:2260
virtual HYPRE_PtrToParSolverFcn SolveFcn() const =0
hypre&#39;s internal Solve function
void SetRelTol(double rel_tol)
Definition: hypre.cpp:3971
void SetKDim(int dim)
Definition: hypre.cpp:2491
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:2335
double b
Definition: lissajous.cpp:42
int to_int(const std::string &str)
Definition: text.hpp:70
Abort on hypre errors (default in base class)
Definition: hypre.hpp:692
void SetPrintLevel(int print_level)
Definition: hypre.hpp:992
void GetBlocks(Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
Definition: hypre.cpp:981
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
Definition: hypre.cpp:2192
HypreParaSails(MPI_Comm comm)
Definition: hypre.cpp:2606
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
PCG Solve function.
Definition: hypre.hpp:787
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:3236
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.cpp:4077
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:1003
HYPRE_Int Randomize(HYPRE_Int seed)
Set random values.
Definition: hypre.cpp:211
Ignore hypre errors (see e.g. HypreADS)
Definition: hypre.hpp:690
double relax_weight
Damping coefficient (usually &lt;= 1)
Definition: hypre.hpp:599
Parallel smoothers in hypre.
Definition: hypre.hpp:581
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:2989
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:939
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:998
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:850
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:585
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:870
Dynamic 2D array using row-major layout.
Definition: array.hpp:316
virtual void SetOperator(const Operator &op)
Definition: hypre.cpp:1942
virtual HYPRE_PtrToParSolverFcn SetupFcn() const =0
hypre&#39;s internal Setup function
void SetLogging(int logging)
Definition: hypre.cpp:2496
virtual ~HypreADS()
Definition: hypre.cpp:3522
HypreParMatrix()
An empty matrix to be used as a reference to an existing matrix.
Definition: hypre.cpp:298
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:868
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2310
HypreParMatrix * GetData()
Definition: hypre.hpp:873
HYPRE_Int * RowPart()
Returns the row partitioning.
Definition: hypre.hpp:372
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
Definition: hypre.cpp:1900
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:73
HypreParMatrix * EliminateCols(const Array< int > &cols)
Definition: hypre.cpp:1413
HypreParVector * X
Definition: hypre.hpp:700
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1084
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:229
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:3751
PCG solver in hypre.
Definition: hypre.hpp:743
GMRES solver in hypre.
Definition: hypre.hpp:798
signed char OwnsOffd() const
Get offd ownership flag.
Definition: hypre.hpp:354
const HYPRE_Int * RowPart() const
Returns the row partitioning (const version)
Definition: hypre.hpp:380
signed char OwnsColMap() const
Get colmap ownership flag.
Definition: hypre.hpp:356
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3267
HypreParVector * X
Definition: hypre.hpp:587
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:4024
void SetNumModes(int num_eigs)
Definition: hypre.cpp:3957
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:905
HypreParMatrix(hypre_ParCSRMatrix *a, bool owner=true)
Converts hypre&#39;s format to HypreParMatrix.
Definition: hypre.hpp:239
virtual ~HypreDiagScale()
Definition: hypre.hpp:874
int GetNumRows() const
Returns the number of rows in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:407
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.hpp:721
HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A + B.
Definition: hypre.cpp:1627
double max_eig_est
Maximal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:619
A class to initialize the size of a Tensor.
Definition: dtensor.hpp:54
void SetRelTol(double rel_tol)
Definition: hypre.cpp:3655
HypreAME(MPI_Comm comm)
Definition: hypre.cpp:3915
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:903
virtual ~HypreIdentity()
Definition: hypre.hpp:855
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3981
void SetFIRCoefficients(double max_eig)
Compute window and Chebyshev coefficients for given polynomial order.
Definition: hypre.cpp:2019
void SetRandomSeed(int s)
Definition: hypre.hpp:1210
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:1425
signed char OwnsDiag() const
Get diag ownership flag.
Definition: hypre.hpp:352
double * fir_coeffs
Combined coefficients for windowing and Chebyshev polynomials.
Definition: hypre.hpp:626
void Threshold(double threshold=0.0)
Remove values smaller in absolute value than some threshold.
Definition: hypre.cpp:1314
double InnerProduct(HypreParVector *x, HypreParVector *y)
Definition: hypre.cpp:246
HYPRE_Int GlobalSize()
Returns the global number of rows.
Definition: hypre.hpp:125
virtual ~HypreSmoother()
Definition: hypre.cpp:2156
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:1217
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2890
void SetData(double *_data)
Sets the data of the Vector and the hypre_ParVector to _data.
Definition: hypre.cpp:205
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
Definition: hypre.cpp:1235
void SetTaubinOptions(double lambda, double mu, int iter)
Set parameters for Taubin&#39;s lambda-mu method.
Definition: hypre.cpp:1912
HypreParVector * B
Right-hand side and solution vectors.
Definition: hypre.hpp:587
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2744
int dim
Definition: ex24.cpp:43
void SetMassMatrix(HypreParMatrix &M)
Definition: hypre.cpp:4017
HYPRE_Int M() const
Returns the global number of rows.
Definition: hypre.hpp:386
void SetZeroInintialIterate()
non-hypre setting
Definition: hypre.hpp:825
virtual ~HyprePCG()
Definition: hypre.cpp:2424
int GetNumCols() const
Returns the number of columns in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:414
void SetOperator(Operator &A)
Definition: hypre.cpp:3695
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2325
void SetPolyOptions(int poly_order, double poly_fraction)
Set parameters for polynomial smoothing.
Definition: hypre.cpp:1906
void GetNumIterations(int &num_iterations)
Definition: hypre.hpp:773
void CopyColStarts()
Definition: hypre.cpp:898
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:944
int relax_times
Number of relaxation sweeps.
Definition: hypre.hpp:597
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
Definition: hypre.cpp:1276
void SetPrintLevel(int logging)
Definition: hypre.cpp:3987
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:257
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:684
ErrorMode error_mode
How to treat hypre errors.
Definition: hypre.hpp:706
HypreParVector & operator=(double d)
Set constant values.
Definition: hypre.cpp:186
const double alpha
Definition: ex15.cpp:336
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:1477
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:847
void SetPrecondUsageMode(int pcg_mode)
Definition: hypre.cpp:3680
virtual ~HypreParMatrix()
Calls hypre&#39;s destroy function.
Definition: hypre.hpp:547
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:697
ID for class HypreParMatrix.
Definition: operator.hpp:233
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2486
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:975
HYPRE_Int * GetColStarts() const
Return the parallel column partitioning array.
Definition: hypre.hpp:436
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:4066
double omega
SOR parameter (usually in (0,2))
Definition: hypre.hpp:601
HYPRE_Int * GetRowStarts() const
Return the parallel row partitioning array.
Definition: hypre.hpp:431
Base class for solvers.
Definition: operator.hpp:500
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1044
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:784
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:187
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:942
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1000
~HypreParVector()
Calls hypre&#39;s destroy function.
Definition: hypre.cpp:221
void SetNumModes(int num_eigs)
Definition: hypre.hpp:1208
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2348
HypreParVector * Z
Definition: hypre.hpp:589
void SetSingularProblem()
Set this option when solving a curl-curl problem with zero mass term.
Definition: hypre.hpp:1039
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:3808
void Read(MPI_Comm comm, const char *fname)
Reads the matrix from a file.
Definition: hypre.cpp:1441
HypreParVector * V
Temporary vectors.
Definition: hypre.hpp:589
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:1894
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
GMRES Setup function.
Definition: hypre.hpp:831
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:603
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:610
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1086