MFEM  v4.4.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-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_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 #if defined(HYPRE_USING_GPU) && \
36  !(defined(HYPRE_USING_CUDA) || defined(HYPRE_USING_HIP))
37 #error "Unsupported GPU build of HYPRE! Only CUDA and HIP builds are supported."
38 #endif
39 #if defined(HYPRE_USING_CUDA) && !defined(MFEM_USE_CUDA)
40 #error "MFEM_USE_CUDA=YES is required when HYPRE is built with CUDA!"
41 #endif
42 #if defined(HYPRE_USING_HIP) && !defined(MFEM_USE_HIP)
43 #error "MFEM_USE_HIP=YES is required when HYPRE is built with HIP!"
44 #endif
45 
46 #include "sparsemat.hpp"
47 #include "hypre_parcsr.hpp"
48 
49 namespace mfem
50 {
51 
52 class ParFiniteElementSpace;
53 class HypreParMatrix;
54 
55 
56 /// @brief A simple singleton class for hypre's global settings, that 1) calls
57 /// HYPRE_Init() and sets some GPU-relevant options at construction and 2) calls
58 /// HYPRE_Finalize() at destruction.
59 class Hypre
60 {
61 public:
62  /// @brief Initialize hypre by calling HYPRE_Init() and set default options.
63  /// After calling Hypre::Init(), hypre will be finalized automatically at
64  /// program exit.
65  ///
66  /// Calling HYPRE_Finalize() manually is not compatible with this class.
67  static void Init() { Instance(); }
68 
69  /// @brief Finalize hypre (called automatically at program exit if
70  /// Hypre::Init() has been called).
71  ///
72  /// Multiple calls to Hypre::Finalize() have no effect. This function can be
73  /// called manually to more precisely control when hypre is finalized.
74  static void Finalize();
75 
76 private:
77  /// Calls HYPRE_Init() when the singleton is constructed.
78  Hypre();
79 
80  /// The singleton destructor (called at program exit) finalizes hypre.
81  ~Hypre() { Finalize(); }
82 
83  /// Set the default hypre global options (mostly GPU-relevant).
84  void SetDefaultOptions();
85 
86  /// Create and return the Hypre singleton object.
87  static Hypre &Instance()
88  {
89  static Hypre hypre;
90  return hypre;
91  }
92 
93  bool finalized = false; ///< Has Hypre::Finalize() been called already?
94 };
95 
96 
97 namespace internal
98 {
99 
100 template <typename int_type>
101 inline int to_int(int_type i)
102 {
103  MFEM_ASSERT(int_type(int(i)) == i, "overflow converting int_type to int");
104  return int(i);
105 }
106 
107 // Specialization for to_int(int)
108 template <> inline int to_int(int i) { return i; }
109 
110 // Convert a HYPRE_Int to int
111 #ifdef HYPRE_BIGINT
112 template <>
113 inline int to_int(HYPRE_Int i)
114 {
115  MFEM_ASSERT(HYPRE_Int(int(i)) == i, "overflow converting HYPRE_Int to int");
116  return int(i);
117 }
118 #endif
119 
120 } // namespace internal
121 
122 
123 /// The MemoryClass used by Hypre objects.
125 {
126 #if !defined(HYPRE_USING_GPU)
127  return MemoryClass::HOST;
128 #elif defined(HYPRE_USING_UNIFIED_MEMORY)
129  return MemoryClass::MANAGED;
130 #else
131  return MemoryClass::DEVICE;
132 #endif
133 }
134 
135 /// The MemoryType used by MFEM when allocating arrays for Hypre objects.
137 {
138 #if !defined(HYPRE_USING_GPU)
139  return Device::GetHostMemoryType();
140 #elif defined(HYPRE_USING_UNIFIED_MEMORY)
141  return MemoryType::MANAGED;
142 #else
143  return MemoryType::DEVICE;
144 #endif
145 }
146 
147 /// Wrapper for hypre's parallel vector class
148 class HypreParVector : public Vector
149 {
150 private:
151  int own_ParVector;
152 
153  /// The actual object
154  hypre_ParVector *x;
155 
156  friend class HypreParMatrix;
157 
158  // Set Vector::data and Vector::size from *x
159  inline void _SetDataAndSize_();
160 
161 public:
162 
163  /// Default constructor, no underlying @a hypre_ParVector is created.
165  {
166  own_ParVector = false;
167  x = NULL;
168  }
169 
170  /** @brief Creates vector with given global size and parallel partitioning of
171  the rows/columns given by @a col. */
172  /** @anchor hypre_partitioning_descr
173  The partitioning is defined in one of two ways depending on the
174  configuration of HYPRE:
175  1. If HYPRE_AssumedPartitionCheck() returns true (the default),
176  then col is of length 2 and the local processor owns columns
177  [col[0],col[1]).
178  2. If HYPRE_AssumedPartitionCheck() returns false, then col is of
179  length (number of processors + 1) and processor P owns columns
180  [col[P],col[P+1]) i.e. each processor has a copy of the same col
181  array. */
182  HypreParVector(MPI_Comm comm, HYPRE_BigInt glob_size, HYPRE_BigInt *col);
183  /** @brief Creates vector with given global size, partitioning of the
184  columns, and data. */
185  /** The data must be allocated and destroyed outside. If @a data_ is NULL, a
186  dummy vector without a valid data array will be created. See @ref
187  hypre_partitioning_descr "here" for a description of the @a col array.
188 
189  If @a is_device_ptr is true, the pointer @a data_ is assumed to be
190  allocated in the memory location HYPRE_MEMORY_DEVICE. */
191  HypreParVector(MPI_Comm comm, HYPRE_BigInt glob_size, double *data_,
192  HYPRE_BigInt *col, bool is_device_ptr = false);
193  /// Creates a deep copy of @a y
194  HypreParVector(const HypreParVector &y);
195  /// Move constructor for HypreParVector. "Steals" data from its argument.
197  /// Creates vector compatible with (i.e. in the domain of) A or A^T
198  explicit HypreParVector(const HypreParMatrix &A, int transpose = 0);
199  /// Creates vector wrapping y
200  explicit HypreParVector(HYPRE_ParVector y);
201  /// Create a true dof parallel vector on a given ParFiniteElementSpace
202  explicit HypreParVector(ParFiniteElementSpace *pfes);
203 
204  /// \brief Constructs a @p HypreParVector *compatible* with the calling vector
205  /// - meaning that it will be the same size and have the same partitioning.
207 
208  /// MPI communicator
209  MPI_Comm GetComm() const { return x->comm; }
210 
211  /// Converts hypre's format to HypreParVector
212  void WrapHypreParVector(hypre_ParVector *y, bool owner=true);
213 
214  /// Returns the parallel row/column partitioning
215  /** See @ref hypre_partitioning_descr "here" for a description of the
216  partitioning array. */
217  inline const HYPRE_BigInt *Partitioning() const { return x->partitioning; }
218 
219  /// @brief Returns a non-const pointer to the parallel row/column
220  /// partitioning.
221  /// Deprecated in favor of HypreParVector::Partitioning() const.
222  MFEM_DEPRECATED
223  inline HYPRE_BigInt *Partitioning() { return x->partitioning; }
224 
225  /// Returns the global number of rows
226  inline HYPRE_BigInt GlobalSize() const { return x->global_size; }
227 
228  /// Typecasting to hypre's hypre_ParVector*
229  operator hypre_ParVector*() const { return x; }
230 #ifndef HYPRE_PAR_VECTOR_STRUCT
231  /// Typecasting to hypre's HYPRE_ParVector, a.k.a. void *
232  operator HYPRE_ParVector() const { return (HYPRE_ParVector) x; }
233 #endif
234  /// Changes the ownership of the vector
235  hypre_ParVector *StealParVector() { own_ParVector = 0; return x; }
236 
237  /// Sets ownership of the internal hypre_ParVector
238  void SetOwnership(int own) { own_ParVector = own; }
239 
240  /// Gets ownership of the internal hypre_ParVector
241  int GetOwnership() const { return own_ParVector; }
242 
243  /// Returns the global vector in each processor
244  Vector* GlobalVector() const;
245 
246  /// Set constant values
247  HypreParVector& operator= (double d);
248  /// Define '=' for hypre vectors.
250  /// Move assignment
252 
253  using Vector::Read;
254 
255  /// Sets the data of the Vector and the hypre_ParVector to @a data_.
256  /** Must be used only for HypreParVector%s that do not own the data,
257  e.g. created with the constructor:
258  HypreParVector(MPI_Comm, HYPRE_BigInt, double *, HYPRE_BigInt *). */
259  void SetData(double *data_);
260 
261  /** @brief Prepare the HypreParVector for read access in hypre's device
262  memory space, HYPRE_MEMORY_DEVICE. */
263  void HypreRead() const;
264 
265  /** @brief Prepare the HypreParVector for read and write access in hypre's
266  device memory space, HYPRE_MEMORY_DEVICE. */
267  void HypreReadWrite();
268 
269  /** @brief Prepare the HypreParVector for write access in hypre's device
270  memory space, HYPRE_MEMORY_DEVICE. */
271  void HypreWrite();
272 
273  /** @brief Replace the HypreParVector's data with the given Memory, @a mem,
274  and prepare the vector for read access in hypre's device memory space,
275  HYPRE_MEMORY_DEVICE. */
276  /** This method must be used with HypreParVector%s that do not own the data,
277  e.g. created with the constructor:
278  HypreParVector(MPI_Comm, HYPRE_BigInt, double *, HYPRE_BigInt *).
279 
280  The Memory @a mem must be accessible with the hypre MemoryClass defined
281  by GetHypreMemoryClass(). */
282  void WrapMemoryRead(const Memory<double> &mem);
283 
284  /** @brief Replace the HypreParVector's data with the given Memory, @a mem,
285  and prepare the vector for read and write access in hypre's device memory
286  space, HYPRE_MEMORY_DEVICE. */
287  /** This method must be used with HypreParVector%s that do not own the data,
288  e.g. created with the constructor:
289  HypreParVector(MPI_Comm, HYPRE_BigInt, double *, HYPRE_BigInt *).
290 
291  The Memory @a mem must be accessible with the hypre MemoryClass defined
292  by GetHypreMemoryClass(). */
294 
295  /** @brief Replace the HypreParVector's data with the given Memory, @a mem,
296  and prepare the vector for write access in hypre's device memory space,
297  HYPRE_MEMORY_DEVICE. */
298  /** This method must be used with HypreParVector%s that do not own the data,
299  e.g. created with the constructor:
300  HypreParVector(MPI_Comm, HYPRE_BigInt, double *, HYPRE_BigInt *).
301 
302  The Memory @a mem must be accessible with the hypre MemoryClass defined
303  by GetHypreMemoryClass(). */
304  void WrapMemoryWrite(Memory<double> &mem);
305 
306  /// Set random values
307  HYPRE_Int Randomize(HYPRE_Int seed);
308 
309  /// Prints the locally owned rows in parallel
310  void Print(const char *fname) const;
311 
312  /// Reads a HypreParVector from files saved with HypreParVector::Print
313  void Read(MPI_Comm comm, const char *fname);
314 
315  /// Calls hypre's destroy function
316  ~HypreParVector();
317 
318 #ifdef MFEM_USE_SUNDIALS
319  /// (DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_PARALLEL.
320  /** @deprecated The returned N_Vector must be destroyed by the caller. */
321  MFEM_DEPRECATED virtual N_Vector ToNVector();
322  using Vector::ToNVector;
323 #endif
324 };
325 
326 /// Returns the inner product of x and y
327 double InnerProduct(HypreParVector &x, HypreParVector &y);
328 double InnerProduct(HypreParVector *x, HypreParVector *y);
329 
330 
331 /** @brief Compute the l_p norm of the Vector which is split without overlap
332  across the given communicator. */
333 double ParNormlp(const Vector &vec, double p, MPI_Comm comm);
334 
335 
336 /// Wrapper for hypre's ParCSR matrix class
337 class HypreParMatrix : public Operator
338 {
339 private:
340  /// The actual object
341  hypre_ParCSRMatrix *A;
342 
343  /// Auxiliary vectors for typecasting
344  mutable HypreParVector *X, *Y;
345  /** @brief Auxiliary buffers for the case when the input or output arrays in
346  methods like Mult(double, const Vector &, double, Vector &) need to be
347  deep copied in order to be used by hypre. */
348  mutable Memory<double> auxX, auxY;
349 
350  // Flags indicating ownership of A->diag->{i,j,data}, A->offd->{i,j,data},
351  // and A->col_map_offd.
352  // The possible values for diagOwner are:
353  // -1: no special treatment of A->diag (default)
354  // when hypre is built with CUDA support, A->diag owns the "host"
355  // pointers (according to A->diag->owns_data)
356  // -2: used when hypre is built with CUDA support, A->diag owns the "hypre"
357  // pointers (according to A->diag->owns_data)
358  // 0: prevent hypre from destroying A->diag->{i,j,data}
359  // 1: same as 0, plus own the "host" A->diag->{i,j}
360  // 2: same as 0, plus own the "host" A->diag->data
361  // 3: same as 0, plus own the "host" A->diag->{i,j,data}
362  // The same values and rules apply to offdOwner and A->offd.
363  // The possible values for colMapOwner are:
364  // -1: no special treatment of A->col_map_offd (default)
365  // 0: prevent hypre from destroying A->col_map_offd
366  // 1: same as 0, plus take ownership of A->col_map_offd
367  // All owned arrays are destroyed with 'delete []'.
368  signed char diagOwner, offdOwner, colMapOwner;
369 
370  // Does the object own the pointer A?
371  signed char ParCSROwner;
372 
373  MemoryIJData mem_diag, mem_offd;
374 
375  // Initialize with defaults. Does not initialize inherited members.
376  void Init();
377 
378  // Delete all owned data. Does not perform re-initialization with defaults.
379  void Destroy();
380 
381  void Read(MemoryClass mc) const;
382  void ReadWrite(MemoryClass mc);
383  // The Boolean flags are used in Destroy().
384  void Write(MemoryClass mc, bool set_diag = true, bool set_offd = true);
385 
386  // Copy (shallow/deep, based on HYPRE_BIGINT) the I and J arrays from csr to
387  // hypre_csr. Shallow copy the data. Return the appropriate ownership flag.
388  // The CSR arrays are wrapped in the mem_csr struct which is used to move
389  // these arrays to device, if necessary.
390  static signed char CopyCSR(SparseMatrix *csr,
391  MemoryIJData &mem_csr,
392  hypre_CSRMatrix *hypre_csr,
393  bool mem_owner);
394  // Copy (shallow or deep, based on HYPRE_BIGINT) the I and J arrays from
395  // bool_csr to hypre_csr. Allocate the data array and set it to all ones.
396  // Return the appropriate ownership flag. The CSR arrays are wrapped in the
397  // mem_csr struct which is used to move these arrays to device, if necessary.
398  static signed char CopyBoolCSR(Table *bool_csr,
399  MemoryIJData &mem_csr,
400  hypre_CSRMatrix *hypre_csr);
401 
402  // Wrap the data from h_mat into mem with the given ownership flag.
403  // If the new Memory arrays in mem are not suitable to be accessed via
404  // GetHypreMemoryClass(), then mem will be re-allocated using the memory type
405  // returned by GetHypreMemoryType(), the data will be deep copied, and h_mat
406  // will be updated with the new pointers.
407  static signed char HypreCsrToMem(hypre_CSRMatrix *h_mat, MemoryType h_mat_mt,
408  bool own_ija, MemoryIJData &mem);
409 
410 public:
411  /// An empty matrix to be used as a reference to an existing matrix
412  HypreParMatrix();
413 
414  /// Converts hypre's format to HypreParMatrix
415  /** If @a owner is false, ownership of @a a is not transferred */
416  void WrapHypreParCSRMatrix(hypre_ParCSRMatrix *a, bool owner = true);
417 
418  /// Converts hypre's format to HypreParMatrix
419  /** If @a owner is false, ownership of @a a is not transferred */
420  explicit HypreParMatrix(hypre_ParCSRMatrix *a, bool owner = true)
421  {
422  Init();
423  WrapHypreParCSRMatrix(a, owner);
424  }
425 
426  /// Creates block-diagonal square parallel matrix.
427  /** Diagonal is given by @a diag which must be in CSR format (finalized). The
428  new HypreParMatrix does not take ownership of any of the input arrays.
429  See @ref hypre_partitioning_descr "here" for a description of the row
430  partitioning array @a row_starts.
431 
432  @warning The ordering of the columns in each row in @a *diag may be
433  changed by this constructor to ensure that the first entry in each row is
434  the diagonal one. This is expected by most hypre functions. */
435  HypreParMatrix(MPI_Comm comm, HYPRE_BigInt glob_size,
436  HYPRE_BigInt *row_starts,
437  SparseMatrix *diag); // constructor with 4 arguments, v1
438 
439  /// Creates block-diagonal rectangular parallel matrix.
440  /** Diagonal is given by @a diag which must be in CSR format (finalized). The
441  new HypreParMatrix does not take ownership of any of the input arrays.
442  See @ref hypre_partitioning_descr "here" for a description of the
443  partitioning arrays @a row_starts and @a col_starts. */
444  HypreParMatrix(MPI_Comm comm, HYPRE_BigInt global_num_rows,
445  HYPRE_BigInt global_num_cols, HYPRE_BigInt *row_starts,
446  HYPRE_BigInt *col_starts,
447  SparseMatrix *diag); // constructor with 6 arguments, v1
448 
449  /// Creates general (rectangular) parallel matrix.
450  /** The new HypreParMatrix does not take ownership of any of the input
451  arrays, if @a own_diag_offd is false (default). If @a own_diag_offd is
452  true, ownership of @a diag and @a offd is transferred to the
453  HypreParMatrix.
454 
455  See @ref hypre_partitioning_descr "here" for a description of the
456  partitioning arrays @a row_starts and @a col_starts. */
457  HypreParMatrix(MPI_Comm comm, HYPRE_BigInt global_num_rows,
458  HYPRE_BigInt global_num_cols, HYPRE_BigInt *row_starts,
459  HYPRE_BigInt *col_starts, SparseMatrix *diag,
460  SparseMatrix *offd, HYPRE_BigInt *cmap,
461  bool own_diag_offd = false); // constructor with 8+1 arguments
462 
463  /// Creates general (rectangular) parallel matrix.
464  /** The new HypreParMatrix takes ownership of all input arrays, except
465  @a col_starts and @a row_starts. See @ref hypre_partitioning_descr "here"
466  for a description of the partitioning arrays @a row_starts and @a
467  col_starts.
468 
469  If @a hypre_arrays is false, all arrays (except @a row_starts and
470  @a col_starts) are assumed to be allocated according to the MemoryType
471  returned by Device::GetHostMemoryType(). If @a hypre_arrays is true, then
472  the same arrays are assumed to be allocated by hypre as host arrays. */
473  HypreParMatrix(MPI_Comm comm,
474  HYPRE_BigInt global_num_rows, HYPRE_BigInt global_num_cols,
475  HYPRE_BigInt *row_starts, HYPRE_BigInt *col_starts,
476  HYPRE_Int *diag_i, HYPRE_Int *diag_j, double *diag_data,
477  HYPRE_Int *offd_i, HYPRE_Int *offd_j, double *offd_data,
478  HYPRE_Int offd_num_cols,
479  HYPRE_BigInt *offd_col_map,
480  bool hypre_arrays = false); // constructor with 13+1 arguments
481 
482  /// Creates a parallel matrix from SparseMatrix on processor 0.
483  /** See @ref hypre_partitioning_descr "here" for a description of the
484  partitioning arrays @a row_starts and @a col_starts. */
485  HypreParMatrix(MPI_Comm comm, HYPRE_BigInt *row_starts,
486  HYPRE_BigInt *col_starts,
487  SparseMatrix *a); // constructor with 4 arguments, v2
488 
489  /// Creates boolean block-diagonal rectangular parallel matrix.
490  /** The new HypreParMatrix does not take ownership of any of the input
491  arrays. See @ref hypre_partitioning_descr "here" for a description of the
492  partitioning arrays @a row_starts and @a col_starts. */
493  HypreParMatrix(MPI_Comm comm, HYPRE_BigInt global_num_rows,
494  HYPRE_BigInt global_num_cols, HYPRE_BigInt *row_starts,
495  HYPRE_BigInt *col_starts,
496  Table *diag); // constructor with 6 arguments, v2
497 
498  /// Creates boolean rectangular parallel matrix.
499  /** The new HypreParMatrix takes ownership of the arrays @a i_diag,
500  @a j_diag, @a i_offd, @a j_offd, and @a cmap; does not take ownership of
501  the arrays @a row and @a col. See @ref hypre_partitioning_descr "here"
502  for a description of the partitioning arrays @a row and @a col. */
503  HypreParMatrix(MPI_Comm comm, int id, int np, HYPRE_BigInt *row,
504  HYPRE_BigInt *col,
505  HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd,
506  HYPRE_Int *j_offd, HYPRE_BigInt *cmap,
507  HYPRE_Int cmap_size); // constructor with 11 arguments
508 
509  /** @brief Creates a general parallel matrix from a local CSR matrix on each
510  processor described by the @a I, @a J and @a data arrays. */
511  /** The local matrix should be of size (local) @a nrows by (global)
512  @a glob_ncols. The new parallel matrix contains copies of all input
513  arrays (so they can be deleted). See @ref hypre_partitioning_descr "here"
514  for a description of the partitioning arrays @a rows and @a cols. */
515  HypreParMatrix(MPI_Comm comm, int nrows, HYPRE_BigInt glob_nrows,
516  HYPRE_BigInt glob_ncols, int *I, HYPRE_BigInt *J,
517  double *data, HYPRE_BigInt *rows,
518  HYPRE_BigInt *cols); // constructor with 9 arguments
519 
520  /** @brief Copy constructor for a ParCSR matrix which creates a deep copy of
521  structure and data from @a P. */
522  HypreParMatrix(const HypreParMatrix &P);
523 
524  /// Make this HypreParMatrix a reference to 'master'
525  void MakeRef(const HypreParMatrix &master);
526 
527  /// MPI communicator
528  MPI_Comm GetComm() const { return A->comm; }
529 
530  /// Typecasting to hypre's hypre_ParCSRMatrix*
531  operator hypre_ParCSRMatrix*() const { return A; }
532 #ifndef HYPRE_PAR_CSR_MATRIX_STRUCT
533  /// Typecasting to hypre's HYPRE_ParCSRMatrix, a.k.a. void *
534  operator HYPRE_ParCSRMatrix() { return (HYPRE_ParCSRMatrix) A; }
535 #endif
536  /// Changes the ownership of the matrix
537  hypre_ParCSRMatrix* StealData();
538 
539  /// Explicitly set the three ownership flags, see docs for diagOwner etc.
540  void SetOwnerFlags(signed char diag, signed char offd, signed char colmap);
541 
542  /// Get diag ownership flag
543  signed char OwnsDiag() const { return diagOwner; }
544  /// Get offd ownership flag
545  signed char OwnsOffd() const { return offdOwner; }
546  /// Get colmap ownership flag
547  signed char OwnsColMap() const { return colMapOwner; }
548 
549  /** If the HypreParMatrix does not own the row-starts array, make a copy of
550  it that the HypreParMatrix will own. If the col-starts array is the same
551  as the row-starts array, col-starts is also replaced. */
552  void CopyRowStarts();
553  /** If the HypreParMatrix does not own the col-starts array, make a copy of
554  it that the HypreParMatrix will own. If the row-starts array is the same
555  as the col-starts array, row-starts is also replaced. */
556  void CopyColStarts();
557 
558  /// Returns the global number of nonzeros
559  inline HYPRE_BigInt NNZ() const { return A->num_nonzeros; }
560  /// Returns the row partitioning
561  /** See @ref hypre_partitioning_descr "here" for a description of the
562  partitioning array. */
563  inline HYPRE_BigInt *RowPart() { return A->row_starts; }
564  /// Returns the column partitioning
565  /** See @ref hypre_partitioning_descr "here" for a description of the
566  partitioning array. */
567  inline HYPRE_BigInt *ColPart() { return A->col_starts; }
568  /// Returns the row partitioning (const version)
569  /** See @ref hypre_partitioning_descr "here" for a description of the
570  partitioning array. */
571  inline const HYPRE_BigInt *RowPart() const { return A->row_starts; }
572  /// Returns the column partitioning (const version)
573  /** See @ref hypre_partitioning_descr "here" for a description of the
574  partitioning array. */
575  inline const HYPRE_BigInt *ColPart() const { return A->col_starts; }
576  /// Returns the global number of rows
577  inline HYPRE_BigInt M() const { return A->global_num_rows; }
578  /// Returns the global number of columns
579  inline HYPRE_BigInt N() const { return A->global_num_cols; }
580 
581  /// Get the local diagonal of the matrix.
582  void GetDiag(Vector &diag) const;
583  /// Get the local diagonal block. NOTE: 'diag' will not own any data.
584  void GetDiag(SparseMatrix &diag) const;
585  /// Get the local off-diagonal block. NOTE: 'offd' will not own any data.
586  void GetOffd(SparseMatrix &offd, HYPRE_BigInt* &cmap) const;
587  /** @brief Get a single SparseMatrix containing all rows from this processor,
588  merged from the diagonal and off-diagonal blocks stored by the
589  HypreParMatrix. */
590  /** @note The number of columns in the SparseMatrix will be the global number
591  of columns in the parallel matrix, so using this method may result in an
592  integer overflow in the column indices. */
593  void MergeDiagAndOffd(SparseMatrix &merged);
594 
595  /// Return the diagonal of the matrix (Operator interface).
596  virtual void AssembleDiagonal(Vector &diag) const { GetDiag(diag); }
597 
598  /** Split the matrix into M x N equally sized blocks of parallel matrices.
599  The size of 'blocks' must already be set to M x N. */
600  void GetBlocks(Array2D<HypreParMatrix*> &blocks,
601  bool interleaved_rows = false,
602  bool interleaved_cols = false) const;
603 
604  /// Returns the transpose of *this
605  HypreParMatrix * Transpose() const;
606 
607  /** Returns principle submatrix given by array of indices of connections
608  with relative size > @a threshold in *this. */
609 #if MFEM_HYPRE_VERSION >= 21800
611  double threshold=0.0) const;
612 #endif
613 
614  /// Returns the number of rows in the diagonal block of the ParCSRMatrix
615  int GetNumRows() const
616  {
617  return internal::to_int(
618  hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(A)));
619  }
620 
621  /// Returns the number of columns in the diagonal block of the ParCSRMatrix
622  int GetNumCols() const
623  {
624  return internal::to_int(
625  hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(A)));
626  }
627 
628  /// Return the global number of rows
630  { return hypre_ParCSRMatrixGlobalNumRows(A); }
631 
632  /// Return the global number of columns
634  { return hypre_ParCSRMatrixGlobalNumCols(A); }
635 
636  /// Return the parallel row partitioning array.
637  /** See @ref hypre_partitioning_descr "here" for a description of the
638  partitioning array. */
639  HYPRE_BigInt *GetRowStarts() const { return hypre_ParCSRMatrixRowStarts(A); }
640 
641  /// Return the parallel column partitioning array.
642  /** See @ref hypre_partitioning_descr "here" for a description of the
643  partitioning array. */
644  HYPRE_BigInt *GetColStarts() const { return hypre_ParCSRMatrixColStarts(A); }
645 
646  virtual MemoryClass GetMemoryClass() const { return GetHypreMemoryClass(); }
647 
648  /// Computes y = alpha * A * x + beta * y
649  HYPRE_Int Mult(HypreParVector &x, HypreParVector &y,
650  double alpha = 1.0, double beta = 0.0) const;
651  /// Computes y = alpha * A * x + beta * y
652  HYPRE_Int Mult(HYPRE_ParVector x, HYPRE_ParVector y,
653  double alpha = 1.0, double beta = 0.0) const;
654  /// Computes y = alpha * A^t * x + beta * y
656  double alpha = 1.0, double beta = 0.0) const;
657 
658  void Mult(double a, const Vector &x, double b, Vector &y) const;
659  void MultTranspose(double a, const Vector &x, double b, Vector &y) const;
660 
661  virtual void Mult(const Vector &x, Vector &y) const
662  { Mult(1.0, x, 0.0, y); }
663  virtual void MultTranspose(const Vector &x, Vector &y) const
664  { MultTranspose(1.0, x, 0.0, y); }
665 
666  /** @brief Computes y = a * |A| * x + b * y, using entry-wise absolute values
667  of the matrix A. */
668  void AbsMult(double a, const Vector &x, double b, Vector &y) const;
669 
670  /** @brief Computes y = a * |At| * x + b * y, using entry-wise absolute
671  values of the transpose of the matrix A. */
672  void AbsMultTranspose(double a, const Vector &x, double b, Vector &y) const;
673 
674  /** @brief The "Boolean" analog of y = alpha * A * x + beta * y, where
675  elements in the sparsity pattern of the matrix are treated as "true". */
676  void BooleanMult(int alpha, const int *x, int beta, int *y)
677  {
678  HostRead();
679  internal::hypre_ParCSRMatrixBooleanMatvec(A, alpha, const_cast<int*>(x),
680  beta, y);
681  HypreRead();
682  }
683 
684  /** @brief The "Boolean" analog of y = alpha * A^T * x + beta * y, where
685  elements in the sparsity pattern of the matrix are treated as "true". */
686  void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
687  {
688  HostRead();
689  internal::hypre_ParCSRMatrixBooleanMatvecT(A, alpha, const_cast<int*>(x),
690  beta, y);
691  HypreRead();
692  }
693 
694  /// Initialize all entries with value.
695  HypreParMatrix &operator=(double value)
696  {
697 #if MFEM_HYPRE_VERSION < 22200
698  internal::hypre_ParCSRMatrixSetConstantValues(A, value);
699 #else
700  hypre_ParCSRMatrixSetConstantValues(A, value);
701 #endif
702  return *this;
703  }
704 
705  /** Perform the operation `*this += B`, assuming that both matrices use the
706  same row and column partitions and the same col_map_offd arrays, or B has
707  an empty off-diagonal block. We also assume that the sparsity pattern of
708  `*this` contains that of `B`. */
709  HypreParMatrix &operator+=(const HypreParMatrix &B) { return Add(1.0, B); }
710 
711  /** Perform the operation `*this += beta*B`, assuming that both matrices use
712  the same row and column partitions and the same col_map_offd arrays, or
713  B has an empty off-diagonal block. We also assume that the sparsity
714  pattern of `*this` contains that of `B`. For a more general case consider
715  the stand-alone function ParAdd described below. */
716  HypreParMatrix &Add(const double beta, const HypreParMatrix &B)
717  {
718  MFEM_VERIFY(internal::hypre_ParCSRMatrixSum(A, beta, B.A) == 0,
719  "error in hypre_ParCSRMatrixSum");
720  return *this;
721  }
722 
723  /** @brief Multiply the HypreParMatrix on the left by a block-diagonal
724  parallel matrix @a D and return the result as a new HypreParMatrix. */
725  /** If @a D has a different number of rows than @a A (this matrix), @a D's
726  row starts array needs to be given (as returned by the methods
727  GetDofOffsets/GetTrueDofOffsets of ParFiniteElementSpace). The new
728  matrix @a D*A uses copies of the row- and column-starts arrays, so "this"
729  matrix and @a row_starts can be deleted.
730  @note This operation is local and does not require communication. */
732  HYPRE_BigInt* row_starts = NULL) const;
733 
734  /// Scale the local row i by s(i).
735  void ScaleRows(const Vector & s);
736  /// Scale the local row i by 1./s(i)
737  void InvScaleRows(const Vector & s);
738  /// Scale all entries by s: A_scaled = s*A.
739  void operator*=(double s);
740 
741  /// Remove values smaller in absolute value than some threshold
742  void Threshold(double threshold = 0.0);
743 
744  /** @brief Wrapper for hypre_ParCSRMatrixDropSmallEntries in different
745  versions of hypre. Drop off-diagonal entries that are smaller than
746  tol * l2 norm of its row */
747  /** For HYPRE versions < 2.14, this method just calls Threshold() with
748  threshold = tol * max(l2 row norm). */
749  void DropSmallEntries(double tol);
750 
751  /// If a row contains only zeros, set its diagonal to 1.
752  void EliminateZeroRows() { hypre_ParCSRMatrixFixZeroRows(A); }
753 
754  /** Eliminate rows and columns from the matrix, and rows from the vector B.
755  Modify B with the BC values in X. */
756  void EliminateRowsCols(const Array<int> &rows_cols, const HypreParVector &X,
757  HypreParVector &B);
758 
759  /** Eliminate rows and columns from the matrix and store the eliminated
760  elements in a new matrix Ae (returned), so that the modified matrix and
761  Ae sum to the original matrix. */
762  HypreParMatrix* EliminateRowsCols(const Array<int> &rows_cols);
763 
764  /** Eliminate columns from the matrix and store the eliminated elements in a
765  new matrix Ae (returned) so that the modified matrix and Ae sum to the
766  original matrix. */
768 
769  /// Eliminate rows from the diagonal and off-diagonal blocks of the matrix.
770  void EliminateRows(const Array<int> &rows);
771 
772  /** @brief Eliminate essential BC specified by @a ess_dof_list from the
773  solution @a X to the r.h.s. @a B. */
774  /** This matrix is the matrix with eliminated BC, while @a Ae is such that
775  (A+Ae) is the original (Neumann) matrix before elimination. */
776  void EliminateBC(const HypreParMatrix &Ae, const Array<int> &ess_dof_list,
777  const Vector &X, Vector &B) const;
778 
779  /// Update the internal hypre_ParCSRMatrix object, A, to be on host.
780  /** After this call A's diagonal and off-diagonal should not be modified
781  until after a suitable call to {Host,Hypre}{Write,ReadWrite}. */
782  void HostRead() const { Read(Device::GetHostMemoryClass()); }
783 
784  /// Update the internal hypre_ParCSRMatrix object, A, to be on host.
785  /** After this call A's diagonal and off-diagonal can be modified on host
786  and subsequent calls to Hypre{Read,Write,ReadWrite} will require a deep
787  copy of the data if hypre is built with device support. */
788  void HostReadWrite() { ReadWrite(Device::GetHostMemoryClass()); }
789 
790  /// Update the internal hypre_ParCSRMatrix object, A, to be on host.
791  /** Similar to HostReadWrite(), except that the data will never be copied
792  from device to host to ensure host contains the correct current data. */
794 
795  /** @brief Update the internal hypre_ParCSRMatrix object, A, to be in hypre
796  memory space. */
797  /** After this call A's diagonal and off-diagonal should not be modified
798  until after a suitable call to {Host,Hypre}{Write,ReadWrite}. */
799  void HypreRead() const { Read(GetHypreMemoryClass()); }
800 
801  /** @brief Update the internal hypre_ParCSRMatrix object, A, to be in hypre
802  memory space. */
803  /** After this call A's diagonal and off-diagonal can be modified in hypre
804  memory space and subsequent calls to Host{Read,Write,ReadWrite} will
805  require a deep copy of the data if hypre is built with device support. */
806  void HypreReadWrite() { ReadWrite(GetHypreMemoryClass()); }
807 
808  /** @brief Update the internal hypre_ParCSRMatrix object, A, to be in hypre
809  memory space. */
810  /** Similar to HostReadWrite(), except that the data will never be copied
811  from host to hypre memory space to ensure the latter contains the correct
812  current data. */
813  void HypreWrite() { Write(GetHypreMemoryClass()); }
814 
815  /// Prints the locally owned rows in parallel
816  void Print(const char *fname, HYPRE_Int offi = 0, HYPRE_Int offj = 0) const;
817  /// Reads the matrix from a file
818  void Read(MPI_Comm comm, const char *fname);
819  /// Read a matrix saved as a HYPRE_IJMatrix
820  void Read_IJMatrix(MPI_Comm comm, const char *fname);
821 
822  /// Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
823  void PrintCommPkg(std::ostream &out = mfem::out) const;
824 
825  /** @brief Print sizes and hashes for all data arrays of the HypreParMatrix
826  from the local MPI rank. */
827  /** This is a compact text representation of the local data of the
828  HypreParMatrix that can be used to compare matrices from different runs
829  without the need to save the whole matrix. */
830  void PrintHash(std::ostream &out) const;
831 
832  /// Calls hypre's destroy function
833  virtual ~HypreParMatrix() { Destroy(); }
834 
835  Type GetType() const { return Hypre_ParCSR; }
836 };
837 
838 #if MFEM_HYPRE_VERSION >= 21800
839 
841 {
842  MATRIX_ONLY,
843  RHS_ONLY,
845 };
846 
847 /** Constructs and applies block diagonal inverse of HypreParMatrix.
848  The enum @a job specifies whether the matrix or the RHS should be
849  scaled (or both). */
850 void BlockInverseScale(const HypreParMatrix *A, HypreParMatrix *C,
851  const Vector *b, HypreParVector *d,
852  int blocksize, BlockInverseScaleJob job);
853 #endif
854 
855 /** @brief Return a new matrix `C = alpha*A + beta*B`, assuming that both `A`
856  and `B` use the same row and column partitions and the same `col_map_offd`
857  arrays. */
858 HypreParMatrix *Add(double alpha, const HypreParMatrix &A,
859  double beta, const HypreParMatrix &B);
860 
861 /** Returns the matrix @a A * @a B. Returned matrix does not necessarily own
862  row or column starts unless the bool @a own_matrix is set to true. */
863 HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B,
864  bool own_matrix = false);
865 /// Returns the matrix A + B
866 /** It is assumed that both matrices use the same row and column partitions and
867  the same col_map_offd arrays. */
868 HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B);
869 
870 /// Returns the matrix P^t * A * P
871 HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P);
872 /// Returns the matrix Rt^t * A * P
873 HypreParMatrix * RAP(const HypreParMatrix * Rt, const HypreParMatrix *A,
874  const HypreParMatrix *P);
875 
876 /// Returns a merged hypre matrix constructed from hypre matrix blocks.
877 /** It is assumed that all block matrices use the same communicator, and the
878  block sizes are consistent in rows and columns. Rows and columns are
879  renumbered but not redistributed in parallel, e.g. the block rows owned by
880  each process remain on that process in the resulting matrix. Some blocks can
881  be NULL. Each block and the entire system can be rectangular. Scalability to
882  extremely large processor counts is limited by global MPI communication, see
883  GatherBlockOffsetData() in hypre.cpp. */
884 HypreParMatrix * HypreParMatrixFromBlocks(Array2D<HypreParMatrix*> &blocks,
885  Array2D<double> *blockCoeff=NULL);
886 
887 /** @brief Eliminate essential BC specified by @a ess_dof_list from the solution
888  @a X to the r.h.s. @a B. */
889 /** Here @a A is a matrix with eliminated BC, while @a Ae is such that (A+Ae) is
890  the original (Neumann) matrix before elimination. */
891 void EliminateBC(const HypreParMatrix &A, const HypreParMatrix &Ae,
892  const Array<int> &ess_dof_list, const Vector &X, Vector &B);
893 
894 
895 /// Parallel smoothers in hypre
896 class HypreSmoother : public Solver
897 {
898 protected:
899  /// The linear system matrix
901  /// Right-hand side and solution vectors
902  mutable HypreParVector *B, *X;
903  /** @brief Auxiliary buffers for the case when the input or output arrays in
904  methods like Mult(const Vector &, Vector &) need to be deep copied in
905  order to be used by hypre. */
907  /// Temporary vectors
908  mutable HypreParVector *V, *Z;
909  /// FIR Filter Temporary Vectors
910  mutable HypreParVector *X0, *X1;
911 
912  /** Smoother type from hypre_ParCSRRelax() in ams.c plus extensions, see the
913  enumeration Type below. */
914  int type;
915  /// Number of relaxation sweeps
917  /// Damping coefficient (usually <= 1)
918  double relax_weight;
919  /// SOR parameter (usually in (0,2))
920  double omega;
921  /// Order of the smoothing polynomial
923  /// Fraction of spectrum to smooth for polynomial relaxation
925  /// Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}
927 
928  /// Taubin's lambda-mu method parameters
929  double lambda;
930  double mu;
932 
933  /// l1 norms of the rows of A
934  double *l1_norms;
935  /// If set, take absolute values of the computed l1_norms
937  /// Number of CG iterations to determine eigenvalue estimates
939  /// Maximal eigenvalue estimate for polynomial smoothing
940  double max_eig_est;
941  /// Minimal eigenvalue estimate for polynomial smoothing
942  double min_eig_est;
943  /// Parameters for windowing function of FIR filter
944  double window_params[3];
945 
946  /// Combined coefficients for windowing and Chebyshev polynomials.
947  double* fir_coeffs;
948 
949  /// A flag that indicates whether the linear system matrix A is symmetric
951 
952 public:
953  /** Hypre smoother types:
954  0 = Jacobi
955  1 = l1-scaled Jacobi
956  2 = l1-scaled block Gauss-Seidel/SSOR
957  4 = truncated l1-scaled block Gauss-Seidel/SSOR
958  5 = lumped Jacobi
959  6 = Gauss-Seidel
960  10 = On-processor forward solve for matrix w/ triangular structure
961  16 = Chebyshev
962  1001 = Taubin polynomial smoother
963  1002 = FIR polynomial smoother. */
964  enum Type { Jacobi = 0, l1Jacobi = 1, l1GS = 2, l1GStr = 4, lumpedJacobi = 5,
965  GS = 6, OPFS = 10, Chebyshev = 16, Taubin = 1001, FIR = 1002
966  };
967 #if !defined(HYPRE_USING_GPU)
968  static constexpr Type default_type = l1GS;
969 #else
970  static constexpr Type default_type = l1Jacobi;
971 #endif
972 
973  HypreSmoother();
974 
976  int relax_times = 1, double relax_weight = 1.0,
977  double omega = 1.0, int poly_order = 2,
978  double poly_fraction = .3, int eig_est_cg_iter = 10);
979 
980  /// Set the relaxation type and number of sweeps
982  /// Set SOR-related parameters
983  void SetSOROptions(double relax_weight, double omega);
984  /// Set parameters for polynomial smoothing
985  /** By default, 10 iterations of CG are used to estimate the eigenvalues.
986  Setting eig_est_cg_iter = 0 uses hypre's hypre_ParCSRMaxEigEstimate() instead. */
987  void SetPolyOptions(int poly_order, double poly_fraction,
988  int eig_est_cg_iter = 10);
989  /// Set parameters for Taubin's lambda-mu method
990  void SetTaubinOptions(double lambda, double mu, int iter);
991 
992  /// Convenience function for setting canonical windowing parameters
993  void SetWindowByName(const char* window_name);
994  /// Set parameters for windowing function for FIR smoother.
995  void SetWindowParameters(double a, double b, double c);
996  /// Compute window and Chebyshev coefficients for given polynomial order.
997  void SetFIRCoefficients(double max_eig);
998 
999  /// After computing l1-norms, replace them with their absolute values.
1000  /** By default, the l1-norms take their sign from the corresponding diagonal
1001  entries in the associated matrix. */
1002  void SetPositiveDiagonal(bool pos = true) { pos_l1_norms = pos; }
1003 
1004  /** Explicitly indicate whether the linear system matrix A is symmetric. If A
1005  is symmetric, the smoother will also be symmetric. In this case, calling
1006  MultTranspose will be redirected to Mult. (This is also done if the
1007  smoother is diagonal.) By default, A is assumed to be nonsymmetric. */
1008  void SetOperatorSymmetry(bool is_sym) { A_is_symmetric = is_sym; }
1009 
1010  /** Set/update the associated operator. Must be called after setting the
1011  HypreSmoother type and options. */
1012  virtual void SetOperator(const Operator &op);
1013 
1014  /// Relax the linear system Ax=b
1015  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
1016  virtual void Mult(const Vector &b, Vector &x) const;
1017 
1018  /// Apply transpose of the smoother to relax the linear system Ax=b
1019  virtual void MultTranspose(const Vector &b, Vector &x) const;
1020 
1021  virtual ~HypreSmoother();
1022 };
1023 
1024 
1025 /// Abstract class for hypre's solvers and preconditioners
1026 class HypreSolver : public Solver
1027 {
1028 public:
1029  /// How to treat errors returned by hypre function calls.
1031  {
1032  IGNORE_HYPRE_ERRORS, ///< Ignore hypre errors (see e.g. HypreADS)
1033  WARN_HYPRE_ERRORS, ///< Issue warnings on hypre errors
1034  ABORT_HYPRE_ERRORS ///< Abort on hypre errors (default in base class)
1035  };
1036 
1037 protected:
1038  /// The linear system matrix
1040 
1041  /// Right-hand side and solution vector
1042  mutable HypreParVector *B, *X;
1043 
1045 
1046  /// Was hypre's Setup function called already?
1047  mutable int setup_called;
1048 
1049  /// How to treat hypre errors.
1051 
1052 public:
1053  HypreSolver();
1054 
1055  HypreSolver(const HypreParMatrix *A_);
1056 
1057  /// Typecast to HYPRE_Solver -- return the solver
1058  virtual operator HYPRE_Solver() const = 0;
1059 
1060  /// hypre's internal Setup function
1061  virtual HYPRE_PtrToParSolverFcn SetupFcn() const = 0;
1062  /// hypre's internal Solve function
1063  virtual HYPRE_PtrToParSolverFcn SolveFcn() const = 0;
1064 
1065  virtual void SetOperator(const Operator &op)
1066  { mfem_error("HypreSolvers do not support SetOperator!"); }
1067 
1068  virtual MemoryClass GetMemoryClass() const { return GetHypreMemoryClass(); }
1069 
1070  /// Solve the linear system Ax=b
1071  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
1072  virtual void Mult(const Vector &b, Vector &x) const;
1073 
1074  /** @brief Set the behavior for treating hypre errors, see the ErrorMode
1075  enum. The default mode in the base class is ABORT_HYPRE_ERRORS. */
1076  /** Currently, there are three cases in derived classes where the error flag
1077  is set to IGNORE_HYPRE_ERRORS:
1078  * in the method HypreBoomerAMG::SetElasticityOptions(), and
1079  * in the constructor of classes HypreAMS and HypreADS.
1080  The reason for this is that a nonzero hypre error is returned) when
1081  hypre_ParCSRComputeL1Norms() encounters zero row in a matrix, which is
1082  expected in some cases with the above solvers. */
1083  void SetErrorMode(ErrorMode err_mode) const { error_mode = err_mode; }
1084 
1085  virtual ~HypreSolver();
1086 };
1087 
1088 
1089 #if MFEM_HYPRE_VERSION >= 21800
1090 /** Preconditioner for HypreParMatrices that are triangular in some ordering.
1091  Finds correct ordering and performs forward substitution on processor
1092  as approximate inverse. Exact on one processor. */
1094 {
1095 public:
1097  explicit HypreTriSolve(const HypreParMatrix &A) : HypreSolver(&A) { }
1098  virtual operator HYPRE_Solver() const { return NULL; }
1099 
1100  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1101  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSROnProcTriSetup; }
1102  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1103  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSROnProcTriSolve; }
1104 
1105  const HypreParMatrix* GetData() const { return A; }
1106 
1107  /// Deprecated. Use HypreTriSolve::GetData() const instead.
1108  MFEM_DEPRECATED HypreParMatrix* GetData()
1109  { return const_cast<HypreParMatrix*>(A); }
1110 
1111  virtual ~HypreTriSolve() { }
1112 };
1113 #endif
1114 
1115 /// PCG solver in hypre
1116 class HyprePCG : public HypreSolver
1117 {
1118 private:
1119  HYPRE_Solver pcg_solver;
1120 
1121  HypreSolver * precond;
1122 
1123 public:
1124  HyprePCG(MPI_Comm comm);
1125 
1126  HyprePCG(const HypreParMatrix &A_);
1127 
1128  virtual void SetOperator(const Operator &op);
1129 
1130  void SetTol(double tol);
1131  void SetAbsTol(double atol);
1132  void SetMaxIter(int max_iter);
1133  void SetLogging(int logging);
1134  void SetPrintLevel(int print_lvl);
1135 
1136  /// Set the hypre solver to be used as a preconditioner
1137  void SetPreconditioner(HypreSolver &precond);
1138 
1139  /** Use the L2 norm of the residual for measuring PCG convergence, plus
1140  (optionally) 1) periodically recompute true residuals from scratch; and
1141  2) enable residual-based stopping criteria. */
1142  void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0);
1143 
1144  /// deprecated: use SetZeroInitialIterate()
1145  MFEM_DEPRECATED void SetZeroInintialIterate() { iterative_mode = false; }
1146 
1147  /// non-hypre setting
1149 
1150  void GetNumIterations(int &num_iterations) const
1151  {
1152  HYPRE_Int num_it;
1153  HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_it);
1154  num_iterations = internal::to_int(num_it);
1155  }
1156 
1157  /// The typecast to HYPRE_Solver returns the internal pcg_solver
1158  virtual operator HYPRE_Solver() const { return pcg_solver; }
1159 
1160  /// PCG Setup function
1161  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1162  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSetup; }
1163  /// PCG Solve function
1164  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1165  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSolve; }
1166 
1167  /// Solve Ax=b with hypre's PCG
1168  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
1169  using HypreSolver::Mult;
1170 
1171  virtual ~HyprePCG();
1172 };
1173 
1174 /// GMRES solver in hypre
1175 class HypreGMRES : public HypreSolver
1176 {
1177 private:
1178  HYPRE_Solver gmres_solver;
1179 
1180  HypreSolver * precond;
1181 
1182  /// Default, generally robust, GMRES options
1183  void SetDefaultOptions();
1184 
1185 public:
1186  HypreGMRES(MPI_Comm comm);
1187 
1188  HypreGMRES(const HypreParMatrix &A_);
1189 
1190  virtual void SetOperator(const Operator &op);
1191 
1192  void SetTol(double tol);
1193  void SetAbsTol(double tol);
1194  void SetMaxIter(int max_iter);
1195  void SetKDim(int dim);
1196  void SetLogging(int logging);
1197  void SetPrintLevel(int print_lvl);
1198 
1199  /// Set the hypre solver to be used as a preconditioner
1200  void SetPreconditioner(HypreSolver &precond);
1201 
1202  /// deprecated: use SetZeroInitialIterate()
1203  MFEM_DEPRECATED void SetZeroInintialIterate() { iterative_mode = false; }
1204 
1205  /// non-hypre setting
1207 
1208  /// The typecast to HYPRE_Solver returns the internal gmres_solver
1209  virtual operator HYPRE_Solver() const { return gmres_solver; }
1210 
1211  /// GMRES Setup function
1212  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1213  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSetup; }
1214  /// GMRES Solve function
1215  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1216  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSolve; }
1217 
1218  /// Solve Ax=b with hypre's GMRES
1219  virtual void Mult (const HypreParVector &b, HypreParVector &x) const;
1220  using HypreSolver::Mult;
1221 
1222  virtual ~HypreGMRES();
1223 };
1224 
1225 /// Flexible GMRES solver in hypre
1226 class HypreFGMRES : public HypreSolver
1227 {
1228 private:
1229  HYPRE_Solver fgmres_solver;
1230 
1231  HypreSolver * precond;
1232 
1233  /// Default, generally robust, FGMRES options
1234  void SetDefaultOptions();
1235 
1236 public:
1237  HypreFGMRES(MPI_Comm comm);
1238 
1239  HypreFGMRES(const HypreParMatrix &A_);
1240 
1241  virtual void SetOperator(const Operator &op);
1242 
1243  void SetTol(double tol);
1244  void SetMaxIter(int max_iter);
1245  void SetKDim(int dim);
1246  void SetLogging(int logging);
1247  void SetPrintLevel(int print_lvl);
1248 
1249  /// Set the hypre solver to be used as a preconditioner
1250  void SetPreconditioner(HypreSolver &precond);
1251 
1252  /// deprecated: use SetZeroInitialIterate()
1253  MFEM_DEPRECATED void SetZeroInintialIterate() { iterative_mode = false; }
1254 
1255  /// non-hypre setting
1257 
1258  /// The typecast to HYPRE_Solver returns the internal fgmres_solver
1259  virtual operator HYPRE_Solver() const { return fgmres_solver; }
1260 
1261  /// FGMRES Setup function
1262  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1263  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRFlexGMRESSetup; }
1264  /// FGMRES Solve function
1265  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1266  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRFlexGMRESSolve; }
1267 
1268  /// Solve Ax=b with hypre's FGMRES
1269  virtual void Mult (const HypreParVector &b, HypreParVector &x) const;
1270  using HypreSolver::Mult;
1271 
1272  virtual ~HypreFGMRES();
1273 };
1274 
1275 /// The identity operator as a hypre solver
1277 {
1278 public:
1279  virtual operator HYPRE_Solver() const { return NULL; }
1280 
1281  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1282  { return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentitySetup; }
1283  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1284  { return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentity; }
1285 
1286  virtual ~HypreIdentity() { }
1287 };
1288 
1289 /// Jacobi preconditioner in hypre
1291 {
1292 public:
1294  explicit HypreDiagScale(const HypreParMatrix &A) : HypreSolver(&A) { }
1295  virtual operator HYPRE_Solver() const { return NULL; }
1296 
1297  virtual void SetOperator(const Operator &op);
1298 
1299  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1300  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScaleSetup; }
1301  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1302  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScale; }
1303 
1304  const HypreParMatrix* GetData() const { return A; }
1305 
1306  /// Deprecated. Use HypreDiagScale::GetData() const instead.
1307  MFEM_DEPRECATED HypreParMatrix* GetData()
1308  { return const_cast<HypreParMatrix*>(A); }
1309 
1310  virtual ~HypreDiagScale() { }
1311 };
1312 
1313 /// The ParaSails preconditioner in hypre
1315 {
1316 private:
1317  HYPRE_Solver sai_precond;
1318 
1319  /// Default, generally robust, ParaSails options
1320  void SetDefaultOptions();
1321 
1322  // If sai_precond is NULL, this method allocates it and sets default options.
1323  // Otherwise the method saves the options from sai_precond, destroys it,
1324  // allocates a new object, and sets its options to the saved values.
1325  void ResetSAIPrecond(MPI_Comm comm);
1326 
1327 public:
1328  HypreParaSails(MPI_Comm comm);
1329 
1331 
1332  virtual void SetOperator(const Operator &op);
1333 
1334  void SetSymmetry(int sym);
1335 
1336  /// The typecast to HYPRE_Solver returns the internal sai_precond
1337  virtual operator HYPRE_Solver() const { return sai_precond; }
1338 
1339  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1340  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSetup; }
1341  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1342  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSolve; }
1343 
1344  virtual ~HypreParaSails();
1345 };
1346 
1347 /** The Euclid preconditioner in Hypre
1348 
1349  Euclid implements the Parallel Incomplete LU factorization technique. For
1350  more information see:
1351 
1352  "A Scalable Parallel Algorithm for Incomplete Factor Preconditioning" by
1353  David Hysom and Alex Pothen, https://doi.org/10.1137/S1064827500376193
1354 */
1355 class HypreEuclid : public HypreSolver
1356 {
1357 private:
1358  HYPRE_Solver euc_precond;
1359 
1360  /// Default, generally robust, Euclid options
1361  void SetDefaultOptions();
1362 
1363  // If euc_precond is NULL, this method allocates it and sets default options.
1364  // Otherwise the method saves the options from euc_precond, destroys it,
1365  // allocates a new object, and sets its options to the saved values.
1366  void ResetEuclidPrecond(MPI_Comm comm);
1367 
1368 public:
1369  HypreEuclid(MPI_Comm comm);
1370 
1371  HypreEuclid(const HypreParMatrix &A);
1372 
1373  virtual void SetOperator(const Operator &op);
1374 
1375  /// The typecast to HYPRE_Solver returns the internal euc_precond
1376  virtual operator HYPRE_Solver() const { return euc_precond; }
1377 
1378  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1379  { return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSetup; }
1380  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1381  { return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSolve; }
1382 
1383  virtual ~HypreEuclid();
1384 };
1385 
1386 #if MFEM_HYPRE_VERSION >= 21900
1387 /**
1388 @brief Wrapper for Hypre's native parallel ILU preconditioner.
1389 
1390 The default ILU factorization type is ILU(k). If you need to change this, or
1391 any other option, you can use the HYPRE_Solver method to cast the object for use
1392 with Hypre's native functions. For example, if want to use natural ordering
1393 rather than RCM reordering, you can use the following approach:
1394 
1395 @code
1396 mfem::HypreILU ilu();
1397 int reorder_type = 0;
1398 HYPRE_ILUSetLocalReordering(ilu, reorder_type);
1399 @endcode
1400 */
1401 class HypreILU : public HypreSolver
1402 {
1403 private:
1404  HYPRE_Solver ilu_precond;
1405 
1406  /// Set the ILU default options
1407  void SetDefaultOptions();
1408 
1409  /** Reset the ILU preconditioner.
1410  @note If ilu_precond is NULL, this method allocates; otherwise it destroys
1411  ilu_precond and allocates a new object. In both cases the default options
1412  are set. */
1413  void ResetILUPrecond();
1414 
1415 public:
1416  /// Constructor; sets the default options
1417  HypreILU();
1418 
1419  virtual ~HypreILU();
1420 
1421  /// Set the fill level for ILU(k); the default is k=1.
1422  void SetLevelOfFill(HYPRE_Int lev_fill);
1423 
1424  /// Set the print level: 0 = none, 1 = setup, 2 = solve, 3 = setup+solve
1425  void SetPrintLevel(HYPRE_Int print_level);
1426 
1427  /// The typecast to HYPRE_Solver returns the internal ilu_precond
1428  virtual operator HYPRE_Solver() const { return ilu_precond; }
1429 
1430  virtual void SetOperator(const Operator &op);
1431 
1432  /// ILU Setup function
1433  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1434  { return (HYPRE_PtrToParSolverFcn) HYPRE_ILUSetup; }
1435 
1436  /// ILU Solve function
1437  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1438  { return (HYPRE_PtrToParSolverFcn) HYPRE_ILUSolve; }
1439 };
1440 #endif
1441 
1442 /// The BoomerAMG solver in hypre
1444 {
1445 private:
1446  HYPRE_Solver amg_precond;
1447 
1448  /// Rigid body modes
1450 
1451  /// Finite element space for elasticity problems, see SetElasticityOptions()
1452  ParFiniteElementSpace *fespace;
1453 
1454  /// Recompute the rigid-body modes vectors (in the rbms array)
1455  void RecomputeRBMs();
1456 
1457  /// Default, generally robust, BoomerAMG options
1458  void SetDefaultOptions();
1459 
1460  // If amg_precond is NULL, allocates it and sets default options.
1461  // Otherwise saves the options from amg_precond, destroys it, allocates a new
1462  // one, and sets its options to the saved values.
1463  void ResetAMGPrecond();
1464 
1465 public:
1466  HypreBoomerAMG();
1467 
1469 
1470  virtual void SetOperator(const Operator &op);
1471 
1472  /** More robust options for systems, such as elasticity. */
1473  void SetSystemsOptions(int dim, bool order_bynodes=false);
1474 
1475  /** A special elasticity version of BoomerAMG that takes advantage of
1476  geometric rigid body modes and could perform better on some problems, see
1477  "Improving algebraic multigrid interpolation operators for linear
1478  elasticity problems", Baker, Kolev, Yang, NLAA 2009, DOI:10.1002/nla.688.
1479  This solver assumes Ordering::byVDIM in the FiniteElementSpace used to
1480  construct A. */
1482 
1483 #if MFEM_HYPRE_VERSION >= 21800
1484  /** Hypre parameters to use AIR AMG solve for advection-dominated problems.
1485  See "Nonsymmetric Algebraic Multigrid Based on Local Approximate Ideal
1486  Restriction (AIR)," Manteuffel, Ruge, Southworth, SISC (2018),
1487  DOI:/10.1137/17M1144350. Options: "distanceR" -> distance of neighbor
1488  DOFs for the restriction operator; options include 1, 2, and 15 (1.5).
1489  Strings "prerelax" and "postrelax" indicate points to relax on:
1490  F = F-points, C = C-points, A = all points. E.g., FFC -> relax on
1491  F-points, relax again on F-points, then relax on C-points. */
1492  void SetAdvectiveOptions(int distance=15, const std::string &prerelax="",
1493  const std::string &postrelax="FFC");
1494 
1495  /// Expert option - consult hypre documentation/team
1496  void SetStrongThresholdR(double strengthR)
1497  { HYPRE_BoomerAMGSetStrongThresholdR(amg_precond, strengthR); }
1498 
1499  /// Expert option - consult hypre documentation/team
1500  void SetFilterThresholdR(double filterR)
1501  { HYPRE_BoomerAMGSetFilterThresholdR(amg_precond, filterR); }
1502 
1503  /// Expert option - consult hypre documentation/team
1504  void SetRestriction(int restrict_type)
1505  { HYPRE_BoomerAMGSetRestriction(amg_precond, restrict_type); }
1506 
1507  /// Expert option - consult hypre documentation/team
1509  { HYPRE_BoomerAMGSetIsTriangular(amg_precond, 1); }
1510 
1511  /// Expert option - consult hypre documentation/team
1512  void SetGMRESSwitchR(int gmres_switch)
1513  { HYPRE_BoomerAMGSetGMRESSwitchR(amg_precond, gmres_switch); }
1514 
1515  /// Expert option - consult hypre documentation/team
1516  void SetCycleNumSweeps(int prerelax, int postrelax)
1517  {
1518  HYPRE_BoomerAMGSetCycleNumSweeps(amg_precond, prerelax, 1);
1519  HYPRE_BoomerAMGSetCycleNumSweeps(amg_precond, postrelax, 2);
1520  }
1521 #endif
1522 
1523  void SetPrintLevel(int print_level)
1524  { HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level); }
1525 
1526  void SetMaxIter(int max_iter)
1527  { HYPRE_BoomerAMGSetMaxIter(amg_precond, max_iter); }
1528 
1529  /// Expert option - consult hypre documentation/team
1530  void SetMaxLevels(int max_levels)
1531  { HYPRE_BoomerAMGSetMaxLevels(amg_precond, max_levels); }
1532 
1533  /// Expert option - consult hypre documentation/team
1534  void SetTol(double tol)
1535  { HYPRE_BoomerAMGSetTol(amg_precond, tol); }
1536 
1537  /// Expert option - consult hypre documentation/team
1538  void SetStrengthThresh(double strength)
1539  { HYPRE_BoomerAMGSetStrongThreshold(amg_precond, strength); }
1540 
1541  /// Expert option - consult hypre documentation/team
1542  void SetInterpolation(int interp_type)
1543  { HYPRE_BoomerAMGSetInterpType(amg_precond, interp_type); }
1544 
1545  /// Expert option - consult hypre documentation/team
1546  void SetCoarsening(int coarsen_type)
1547  { HYPRE_BoomerAMGSetCoarsenType(amg_precond, coarsen_type); }
1548 
1549  /// Expert option - consult hypre documentation/team
1550  void SetRelaxType(int relax_type)
1551  { HYPRE_BoomerAMGSetRelaxType(amg_precond, relax_type); }
1552 
1553  /// Expert option - consult hypre documentation/team
1554  void SetCycleType(int cycle_type)
1555  { HYPRE_BoomerAMGSetCycleType(amg_precond, cycle_type); }
1556 
1557  void GetNumIterations(int &num_iterations) const
1558  {
1559  HYPRE_Int num_it;
1560  HYPRE_BoomerAMGGetNumIterations(amg_precond, &num_it);
1561  num_iterations = internal::to_int(num_it);
1562  }
1563 
1564  /// Expert option - consult hypre documentation/team
1565  void SetNodal(int blocksize)
1566  {
1567  HYPRE_BoomerAMGSetNumFunctions(amg_precond, blocksize);
1568  HYPRE_BoomerAMGSetNodal(amg_precond, 1);
1569  }
1570 
1571  /// Expert option - consult hypre documentation/team
1572  void SetAggressiveCoarsening(int num_levels)
1573  { HYPRE_BoomerAMGSetAggNumLevels(amg_precond, num_levels); }
1574 
1575  /// The typecast to HYPRE_Solver returns the internal amg_precond
1576  virtual operator HYPRE_Solver() const { return amg_precond; }
1577 
1578  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1579  { return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSetup; }
1580  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1581  { return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSolve; }
1582 
1583  using HypreSolver::Mult;
1584 
1585  virtual ~HypreBoomerAMG();
1586 };
1587 
1588 /// Compute the discrete gradient matrix between the nodal linear and ND1 spaces
1589 HypreParMatrix* DiscreteGrad(ParFiniteElementSpace *edge_fespace,
1590  ParFiniteElementSpace *vert_fespace);
1591 /// Compute the discrete curl matrix between the ND1 and RT0 spaces
1592 HypreParMatrix* DiscreteCurl(ParFiniteElementSpace *face_fespace,
1593  ParFiniteElementSpace *edge_fespace);
1594 
1595 /// The Auxiliary-space Maxwell Solver in hypre
1596 class HypreAMS : public HypreSolver
1597 {
1598 private:
1599  /// Constuct AMS solver from finite element space
1600  void Init(ParFiniteElementSpace *edge_space);
1601 
1602  HYPRE_Solver ams;
1603 
1604  /// Vertex coordinates
1605  HypreParVector *x, *y, *z;
1606  /// Discrete gradient matrix
1607  HypreParMatrix *G;
1608  /// Nedelec interpolation matrix and its components
1609  HypreParMatrix *Pi, *Pix, *Piy, *Piz;
1610 
1611 public:
1612  HypreAMS(ParFiniteElementSpace *edge_fespace);
1613 
1614  HypreAMS(const HypreParMatrix &A, ParFiniteElementSpace *edge_fespace);
1615 
1616  virtual void SetOperator(const Operator &op);
1617 
1618  void SetPrintLevel(int print_lvl);
1619 
1620  /// Set this option when solving a curl-curl problem with zero mass term
1621  void SetSingularProblem() { HYPRE_AMSSetBetaPoissonMatrix(ams, NULL); }
1622 
1623  /// The typecast to HYPRE_Solver returns the internal ams object
1624  virtual operator HYPRE_Solver() const { return ams; }
1625 
1626  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1627  { return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSetup; }
1628  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1629  { return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSolve; }
1630 
1631  virtual ~HypreAMS();
1632 };
1633 
1634 /// The Auxiliary-space Divergence Solver in hypre
1635 class HypreADS : public HypreSolver
1636 {
1637 private:
1638  /// Constuct ADS solver from finite element space
1639  void Init(ParFiniteElementSpace *face_fespace);
1640 
1641  HYPRE_Solver ads;
1642 
1643  /// Vertex coordinates
1644  HypreParVector *x, *y, *z;
1645  /// Discrete gradient matrix
1646  HypreParMatrix *G;
1647  /// Discrete curl matrix
1648  HypreParMatrix *C;
1649  /// Nedelec interpolation matrix and its components
1650  HypreParMatrix *ND_Pi, *ND_Pix, *ND_Piy, *ND_Piz;
1651  /// Raviart-Thomas interpolation matrix and its components
1652  HypreParMatrix *RT_Pi, *RT_Pix, *RT_Piy, *RT_Piz;
1653 
1654 public:
1655  HypreADS(ParFiniteElementSpace *face_fespace);
1656 
1657  HypreADS(const HypreParMatrix &A, ParFiniteElementSpace *face_fespace);
1658 
1659  virtual void SetOperator(const Operator &op);
1660 
1661  void SetPrintLevel(int print_lvl);
1662 
1663  /// The typecast to HYPRE_Solver returns the internal ads object
1664  virtual operator HYPRE_Solver() const { return ads; }
1665 
1666  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
1667  { return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSetup; }
1668  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
1669  { return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSolve; }
1670 
1671  virtual ~HypreADS();
1672 };
1673 
1674 /** LOBPCG eigenvalue solver in hypre
1675 
1676  The Locally Optimal Block Preconditioned Conjugate Gradient (LOBPCG)
1677  eigenvalue solver is designed to find the lowest eigenmodes of the
1678  generalized eigenvalue problem:
1679  A x = lambda M x
1680  where A is symmetric, potentially indefinite and M is symmetric positive
1681  definite. The eigenvectors are M-orthonormal, meaning that
1682  x^T M x = 1 and x^T M y = 0,
1683  if x and y are distinct eigenvectors. The matrix M is optional and is
1684  assumed to be the identity if left unset.
1685 
1686  The efficiency of LOBPCG relies on the availability of a suitable
1687  preconditioner for the matrix A. The preconditioner is supplied through the
1688  SetPreconditioner() method. It should be noted that the operator used with
1689  the preconditioner need not be A itself.
1690 
1691  For more information regarding LOBPCG see "Block Locally Optimal
1692  Preconditioned Eigenvalue Xolvers (BLOPEX) in Hypre and PETSc" by
1693  A. Knyazev, M. Argentati, I. Lashuk, and E. Ovtchinnikov, SISC, 29(5),
1694  2224-2239, 2007.
1695 */
1697 {
1698 private:
1699  MPI_Comm comm;
1700  int myid;
1701  int numProcs;
1702  int nev; // Number of desired eigenmodes
1703  int seed; // Random seed used for initial vectors
1704 
1705  HYPRE_BigInt glbSize; // Global number of DoFs in the linear system
1706  HYPRE_BigInt * part; // Row partitioning of the linear system
1707 
1708  // Pointer to HYPRE's solver struct
1709  HYPRE_Solver lobpcg_solver;
1710 
1711  // Interface for matrix storage type
1712  mv_InterfaceInterpreter interpreter;
1713 
1714  // Interface for setting up and performing matrix-vector products
1715  HYPRE_MatvecFunctions matvec_fn;
1716 
1717  // Eigenvalues
1718  Array<double> eigenvalues;
1719 
1720  // Forward declaration
1721  class HypreMultiVector;
1722 
1723  // MultiVector to store eigenvectors
1724  HypreMultiVector * multi_vec;
1725 
1726  // Empty vectors used to setup the matrices and preconditioner
1727  HypreParVector * x;
1728 
1729  // An optional operator which projects vectors into a desired subspace
1730  Operator * subSpaceProj;
1731 
1732  /// Internal class to represent a set of eigenvectors
1733  class HypreMultiVector
1734  {
1735  private:
1736  // Pointer to hypre's multi-vector object
1737  mv_MultiVectorPtr mv_ptr;
1738 
1739  // Wrappers for each member of the multivector
1740  HypreParVector ** hpv;
1741 
1742  // Number of vectors in the multivector
1743  int nv;
1744 
1745  public:
1746  HypreMultiVector(int n, HypreParVector & v,
1747  mv_InterfaceInterpreter & interpreter);
1748  ~HypreMultiVector();
1749 
1750  /// Set random values
1751  void Randomize(HYPRE_Int seed);
1752 
1753  /// Extract a single HypreParVector object
1754  HypreParVector & GetVector(unsigned int i);
1755 
1756  /// Transfers ownership of data to returned array of vectors
1757  HypreParVector ** StealVectors();
1758 
1759  operator mv_MultiVectorPtr() const { return mv_ptr; }
1760 
1761  mv_MultiVectorPtr & GetMultiVector() { return mv_ptr; }
1762  };
1763 
1764  static void * OperatorMatvecCreate( void *A, void *x );
1765  static HYPRE_Int OperatorMatvec( void *matvec_data,
1766  HYPRE_Complex alpha,
1767  void *A,
1768  void *x,
1769  HYPRE_Complex beta,
1770  void *y );
1771  static HYPRE_Int OperatorMatvecDestroy( void *matvec_data );
1772 
1773  static HYPRE_Int PrecondSolve(void *solver,
1774  void *A,
1775  void *b,
1776  void *x);
1777  static HYPRE_Int PrecondSetup(void *solver,
1778  void *A,
1779  void *b,
1780  void *x);
1781 
1782 public:
1783  HypreLOBPCG(MPI_Comm comm);
1784  ~HypreLOBPCG();
1785 
1786  void SetTol(double tol);
1787  void SetRelTol(double rel_tol);
1788  void SetMaxIter(int max_iter);
1789  void SetPrintLevel(int logging);
1790  void SetNumModes(int num_eigs) { nev = num_eigs; }
1791  void SetPrecondUsageMode(int pcg_mode);
1792  void SetRandomSeed(int s) { seed = s; }
1793  void SetInitialVectors(int num_vecs, HypreParVector ** vecs);
1794 
1795  // The following four methods support general operators
1796  void SetPreconditioner(Solver & precond);
1797  void SetOperator(Operator & A);
1798  void SetMassMatrix(Operator & M);
1799  void SetSubSpaceProjector(Operator & proj) { subSpaceProj = &proj; }
1800 
1801  /// Solve the eigenproblem
1802  void Solve();
1803 
1804  /// Collect the converged eigenvalues
1805  void GetEigenvalues(Array<double> & eigenvalues) const;
1806 
1807  /// Extract a single eigenvector
1808  const HypreParVector & GetEigenvector(unsigned int i) const;
1809 
1810  /// Transfer ownership of the converged eigenvectors
1811  HypreParVector ** StealEigenvectors() { return multi_vec->StealVectors(); }
1812 };
1813 
1814 /** AME eigenvalue solver in hypre
1815 
1816  The Auxiliary space Maxwell Eigensolver (AME) is designed to find
1817  the lowest eigenmodes of the generalized eigenvalue problem:
1818  Curl Curl x = lambda M x
1819  where the Curl Curl operator is discretized using Nedelec finite element
1820  basis functions. Properties of this discretization are essential to
1821  eliminating the large null space of the Curl Curl operator.
1822 
1823  This eigensolver relies upon the LOBPCG eigensolver internally. It is also
1824  expected that the preconditioner supplied to this method will be the
1825  HypreAMS preconditioner defined above.
1826 
1827  As with LOBPCG, the operator set in the preconditioner need not be the same
1828  as A. This flexibility may be useful in solving eigenproblems which bare a
1829  strong resemblance to the Curl Curl problems for which AME is designed.
1830 
1831  Unlike LOBPCG, this eigensolver requires that the mass matrix be set.
1832  It is possible to circumvent this by passing an identity operator as the
1833  mass matrix but it seems unlikely that this would be useful so it is not the
1834  default behavior.
1835 */
1837 {
1838 private:
1839  int myid;
1840  int numProcs;
1841  int nev; // Number of desired eigenmodes
1842  bool setT;
1843 
1844  // Pointer to HYPRE's AME solver struct
1845  HYPRE_Solver ame_solver;
1846 
1847  // Pointer to HYPRE's AMS solver struct
1848  HypreSolver * ams_precond;
1849 
1850  // Eigenvalues
1851  HYPRE_Real * eigenvalues;
1852 
1853  // MultiVector to store eigenvectors
1854  HYPRE_ParVector * multi_vec;
1855 
1856  // HypreParVector wrappers to contain eigenvectors
1857  mutable HypreParVector ** eigenvectors;
1858 
1859  void createDummyVectors() const;
1860 
1861 public:
1862  HypreAME(MPI_Comm comm);
1863  ~HypreAME();
1864 
1865  void SetTol(double tol);
1866  void SetRelTol(double rel_tol);
1867  void SetMaxIter(int max_iter);
1868  void SetPrintLevel(int logging);
1869  void SetNumModes(int num_eigs);
1870 
1871  // The following four methods support operators of type HypreParMatrix.
1872  void SetPreconditioner(HypreSolver & precond);
1873  void SetOperator(const HypreParMatrix & A);
1874  void SetMassMatrix(const HypreParMatrix & M);
1875 
1876  /// Solve the eigenproblem
1877  void Solve();
1878 
1879  /// Collect the converged eigenvalues
1880  void GetEigenvalues(Array<double> & eigenvalues) const;
1881 
1882  /// Extract a single eigenvector
1883  const HypreParVector & GetEigenvector(unsigned int i) const;
1884 
1885  /// Transfer ownership of the converged eigenvectors
1887 };
1888 
1889 }
1890 
1891 #endif // MFEM_USE_MPI
1892 
1893 #endif
Type GetType() const
Definition: hypre.hpp:835
MemoryType GetHypreMemoryType()
The MemoryType used by MFEM when allocating arrays for Hypre objects.
Definition: hypre.hpp:136
virtual ~HypreBoomerAMG()
Definition: hypre.cpp:4949
void SetAbsTol(double atol)
Definition: hypre.cpp:3780
Memory< double > auxX
Definition: hypre.hpp:906
void SetTol(double tol)
Definition: hypre.cpp:3775
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:3990
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:2240
HypreADS(ParFiniteElementSpace *face_fespace)
Definition: hypre.cpp:5214
const HypreParMatrix * GetData() const
Definition: hypre.hpp:1105
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(), hypre will be finalized automatically at program exit.
Definition: hypre.hpp:67
static constexpr Type default_type
Definition: hypre.hpp:968
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:528
HypreEuclid(MPI_Comm comm)
Definition: hypre.cpp:4361
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:1083
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:4110
Vector * GlobalVector() const
Returns the global vector in each processor.
Definition: hypre.cpp:268
void EliminateBC(const HypreParMatrix &A, const HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
Eliminate essential BC specified by ess_dof_list from the solution X to the r.h.s. B.
Definition: hypre.cpp:3053
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1596
HypreParVector * X0
FIR Filter Temporary Vectors.
Definition: hypre.hpp:910
void SetStrengthThresh(double strength)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1538
double min_eig_est
Minimal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:942
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1283
static MemoryClass GetHostMemoryClass()
Get the current Host MemoryClass. This is the MemoryClass used by most MFEM host Memory objects...
Definition: device.hpp:268
ErrorMode
How to treat errors returned by hypre function calls.
Definition: hypre.hpp:1030
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:1635
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
FGMRES Solve function.
Definition: hypre.hpp:1265
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:5502
MFEM_DEPRECATED HypreParMatrix * GetData()
Deprecated. Use HypreTriSolve::GetData() const instead.
Definition: hypre.hpp:1108
void SetMassMatrix(const HypreParMatrix &M)
Definition: hypre.cpp:5979
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:2673
virtual ~HypreEuclid()
Definition: hypre.cpp:4426
void WrapMemoryWrite(Memory< double > &mem)
Replace the HypreParVector&#39;s data with the given Memory, mem, and prepare the vector for write access...
Definition: hypre.cpp:372
int setup_called
Was hypre&#39;s Setup function called already?
Definition: hypre.hpp:1047
Device memory; using CUDA or HIP *Malloc and *Free.
void HypreRead() const
Prepare the HypreParVector for read access in hypre&#39;s device memory space, HYPRE_MEMORY_DEVICE.
Definition: hypre.cpp:315
void EliminateBC(const HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B) const
Eliminate essential BC specified by ess_dof_list from the solution X to the r.h.s. B.
Definition: hypre.cpp:2292
void HypreRead() const
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
Definition: hypre.hpp:799
HYPRE_BigInt GlobalSize() const
Returns the global number of rows.
Definition: hypre.hpp:226
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:3753
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:4157
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to &#39;master&#39;.
Definition: hypre.cpp:1345
void Read_IJMatrix(MPI_Comm comm, const char *fname)
Read a matrix saved as a HYPRE_IJMatrix.
Definition: hypre.cpp:2373
static void Finalize()
Finalize hypre (called automatically at program exit if Hypre::Init() has been called).
Definition: hypre.cpp:45
HypreParMatrix & Add(const double beta, const HypreParMatrix &B)
Definition: hypre.hpp:716
HypreParVector CreateCompatibleVector() const
Constructs a HypreParVector compatible with the calling vector.
Definition: hypre.cpp:242
HypreParVector * B
Right-hand side and solution vector.
Definition: hypre.hpp:1042
Wrapper for Hypre&#39;s native parallel ILU preconditioner.
Definition: hypre.hpp:1401
void GetNumIterations(int &num_iterations) const
Definition: hypre.hpp:1150
const HypreParMatrix * GetData() const
Definition: hypre.hpp:1304
double window_params[3]
Parameters for windowing function of FIR filter.
Definition: hypre.hpp:944
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: hypre.hpp:661
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1628
void SetInitialVectors(int num_vecs, HypreParVector **vecs)
Definition: hypre.cpp:5728
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:3938
void SetWindowByName(const char *window_name)
Convenience function for setting canonical windowing parameters.
Definition: hypre.cpp:3252
HypreTriSolve(const HypreParMatrix &A)
Definition: hypre.hpp:1097
Issue warnings on hypre errors.
Definition: hypre.hpp:1033
void SetPreconditioner(HypreSolver &precond)
Definition: hypre.cpp:5958
HyprePCG(MPI_Comm comm)
Definition: hypre.cpp:3735
void SetTol(double tol)
Definition: hypre.cpp:5927
HypreILU()
Constructor; sets the default options.
Definition: hypre.cpp:4433
void SetRestriction(int restrict_type)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1504
HypreAMS(ParFiniteElementSpace *edge_fespace)
Definition: hypre.cpp:4959
void SetKDim(int dim)
Definition: hypre.cpp:4142
virtual ~HypreGMRES()
Definition: hypre.cpp:4072
void SetCycleNumSweeps(int prerelax, int postrelax)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1516
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:4152
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:5462
void SetPreconditioner(Solver &precond)
Definition: hypre.cpp:5644
MFEM_DEPRECATED void SetZeroInintialIterate()
deprecated: use SetZeroInitialIterate()
Definition: hypre.hpp:1203
void SetMassMatrix(Operator &M)
Definition: hypre.cpp:5700
Abstract parallel finite element space.
Definition: pfespace.hpp:28
HypreParVector()
Default constructor, no underlying hypre_ParVector is created.
Definition: hypre.hpp:164
void SetPrintLevel(int logging)
Definition: hypre.cpp:5629
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:655
virtual ~HypreAMS()
Definition: hypre.cpp:5194
void SetTol(double tol)
Definition: hypre.cpp:3960
virtual ~HypreILU()
Definition: hypre.cpp:4505
void SetAdvectiveOptions(int distance=15, const std::string &prerelax="", const std::string &postrelax="FFC")
Definition: hypre.cpp:4841
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
GMRES Solve function.
Definition: hypre.hpp:1215
void BlockInverseScale(const HypreParMatrix *A, HypreParMatrix *C, const Vector *b, HypreParVector *d, int blocksize, BlockInverseScaleJob job)
Definition: hypre.cpp:2566
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.hpp:1811
void HostWrite()
Update the internal hypre_ParCSRMatrix object, A, to be on host.
Definition: hypre.hpp:793
HypreParMatrix & operator+=(const HypreParMatrix &B)
Definition: hypre.hpp:709
bool pos_l1_norms
If set, take absolute values of the computed l1_norms.
Definition: hypre.hpp:936
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
Definition: hypre.cpp:2652
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3985
void SetZeroInitialIterate()
non-hypre setting
Definition: hypre.hpp:1256
void SetStrongThresholdR(double strengthR)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1496
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:209
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
Definition: hypre.cpp:1981
virtual MFEM_DEPRECATED N_Vector ToNVector()
(DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_PARALLEL.
Definition: hypre.cpp:418
double poly_fraction
Fraction of spectrum to smooth for polynomial relaxation.
Definition: hypre.hpp:924
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3795
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:4486
HypreParMatrix * DiscreteGrad(ParFiniteElementSpace *edge_fespace, ParFiniteElementSpace *vert_fespace)
Compute the discrete gradient matrix between the nodal linear and ND1 spaces.
int poly_scale
Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}.
Definition: hypre.hpp:926
void MergeDiagAndOffd(SparseMatrix &merged)
Get a single SparseMatrix containing all rows from this processor, merged from the diagonal and off-d...
Definition: hypre.cpp:1559
void GetEigenvalues(Array< double > &eigenvalues) const
Collect the converged eigenvalues.
Definition: hypre.cpp:5710
void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
The &quot;Boolean&quot; analog of y = alpha * A^T * x + beta * y, where elements in the sparsity pattern of the...
Definition: hypre.hpp:686
void SetTol(double tol)
Definition: hypre.cpp:5607
HypreParVector * X1
Definition: hypre.hpp:910
HypreFGMRES(MPI_Comm comm)
Definition: hypre.cpp:4078
HYPRE_BigInt N() const
Returns the global number of columns.
Definition: hypre.hpp:579
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s FGMRES.
Definition: hypre.cpp:4166
HYPRE_BigInt M() const
Returns the global number of rows.
Definition: hypre.hpp:577
The identity operator as a hypre solver.
Definition: hypre.hpp:1276
void SetWindowParameters(double a, double b, double c)
Set parameters for windowing function for FIR smoother.
Definition: hypre.cpp:3267
void EliminateZeroRows()
If a row contains only zeros, set its diagonal to 1.
Definition: hypre.hpp:752
virtual ~HypreTriSolve()
Definition: hypre.hpp:1111
HypreGMRES(MPI_Comm comm)
Definition: hypre.cpp:3906
void GetOffd(SparseMatrix &offd, HYPRE_BigInt *&cmap) const
Get the local off-diagonal block. NOTE: &#39;offd&#39; will not own any data.
Definition: hypre.cpp:1553
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:4326
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s GMRES.
Definition: hypre.cpp:4000
virtual void AssembleDiagonal(Vector &diag) const
Return the diagonal of the matrix (Operator interface).
Definition: hypre.hpp:596
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
Definition: hypre.hpp:1002
void Print(const char *fname) const
Prints the locally owned rows in parallel.
Definition: hypre.cpp:391
void HypreWrite()
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
Definition: hypre.hpp:813
void AbsMultTranspose(double a, const Vector &x, double b, Vector &y) const
Computes y = a * |At| * x + b * y, using entry-wise absolute values of the transpose of the matrix A...
Definition: hypre.cpp:1862
HypreLOBPCG(MPI_Comm comm)
Definition: hypre.cpp:5577
HYPRE_BigInt * GetRowStarts() const
Return the parallel row partitioning array.
Definition: hypre.hpp:639
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:1916
MFEM_DEPRECATED HYPRE_BigInt * Partitioning()
Returns a non-const pointer to the parallel row/column partitioning. Deprecated in favor of HypreParV...
Definition: hypre.hpp:223
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
Definition: hypre.cpp:3426
hypre_ParVector * StealParVector()
Changes the ownership of the vector.
Definition: hypre.hpp:235
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1443
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:4244
void SetRelaxType(int relax_type)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1550
void SetSymmetry(int sym)
Definition: hypre.cpp:4350
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
Definition: hypre.hpp:646
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:663
virtual ~HypreParaSails()
Definition: hypre.cpp:4355
The ParaSails preconditioner in hypre.
Definition: hypre.hpp:1314
void SetMaxIter(int max_iter)
Definition: hypre.cpp:5623
double * l1_norms
l1 norms of the rows of A
Definition: hypre.hpp:934
Data type sparse matrix.
Definition: sparsemat.hpp:46
BlockInverseScaleJob
Definition: hypre.hpp:840
void SetOwnerFlags(signed char diag, signed char offd, signed char colmap)
Explicitly set the three ownership flags, see docs for diagOwner etc.
Definition: hypre.cpp:1382
void PrintHash(std::ostream &out) const
Print sizes and hashes for all data arrays of the HypreParMatrix from the local MPI rank...
Definition: hypre.cpp:2431
void SetLogging(int logging)
Definition: hypre.cpp:3790
Jacobi preconditioner in hypre.
Definition: hypre.hpp:1290
virtual ~HypreSolver()
Definition: hypre.cpp:3726
virtual HYPRE_PtrToParSolverFcn SolveFcn() const =0
hypre&#39;s internal Solve function
void SetRelTol(double rel_tol)
Definition: hypre.cpp:5933
void SetKDim(int dim)
Definition: hypre.cpp:3975
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:154
void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0)
Definition: hypre.cpp:3810
double b
Definition: lissajous.cpp:42
int to_int(const std::string &str)
Convert a string to an int.
Definition: text.hpp:62
Abort on hypre errors (default in base class)
Definition: hypre.hpp:1034
void SetLogging(int logging)
Definition: hypre.cpp:4147
void SetPrintLevel(int print_level)
Definition: hypre.hpp:1523
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
FGMRES Setup function.
Definition: hypre.hpp:1262
void SetOperatorSymmetry(bool is_sym)
Definition: hypre.hpp:1008
void GetBlocks(Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
Definition: hypre.cpp:1583
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
Definition: hypre.cpp:3618
HypreParaSails(MPI_Comm comm)
Definition: hypre.cpp:4262
void HypreWrite()
Prepare the HypreParVector for write access in hypre&#39;s device memory space, HYPRE_MEMORY_DEVICE.
Definition: hypre.cpp:334
void SetTol(double tol)
Definition: hypre.cpp:4132
void WrapHypreParVector(hypre_ParVector *y, bool owner=true)
Converts hypre&#39;s format to HypreParVector.
Definition: hypre.cpp:259
void SetTol(double tol)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1534
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
PCG Solve function.
Definition: hypre.hpp:1164
void Read(MPI_Comm comm, const char *fname)
Reads a HypreParVector from files saved with HypreParVector::Print.
Definition: hypre.cpp:396
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:5176
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
Definition: hypre.hpp:1068
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.cpp:6033
const HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:1039
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:1607
HYPRE_Int Randomize(HYPRE_Int seed)
Set random values.
Definition: hypre.cpp:386
void SetFilterThresholdR(double filterR)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1500
Ignore hypre errors (see e.g. HypreADS)
Definition: hypre.hpp:1032
double relax_weight
Damping coefficient (usually &lt;= 1)
Definition: hypre.hpp:918
Parallel smoothers in hypre.
Definition: hypre.hpp:896
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:4793
Memory< double > auxB
Definition: hypre.hpp:1044
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:1480
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1578
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
ILU Solve function.
Definition: hypre.hpp:1437
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1281
void SetZeroInitialIterate()
non-hypre setting
Definition: hypre.hpp:1206
void SetLevelOfFill(HYPRE_Int lev_fill)
Set the fill level for ILU(k); the default is k=1.
Definition: hypre.cpp:4476
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:900
void SetCoarsening(int coarsen_type)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1546
void SetGMRESSwitchR(int gmres_switch)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1512
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1301
Dynamic 2D array using row-major layout.
Definition: array.hpp:351
virtual void SetOperator(const Operator &op)
Definition: hypre.cpp:3274
virtual HYPRE_PtrToParSolverFcn SetupFcn() const =0
hypre&#39;s internal Setup function
void SetLogging(int logging)
Definition: hypre.cpp:3980
virtual ~HypreADS()
Definition: hypre.cpp:5480
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
Definition: hypre.hpp:633
HypreParMatrix()
An empty matrix to be used as a reference to an existing matrix.
Definition: hypre.cpp:614
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1299
bool A_is_symmetric
A flag that indicates whether the linear system matrix A is symmetric.
Definition: hypre.hpp:950
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3785
void SetMaxIter(int max_iter)
Definition: hypre.hpp:1526
void SetMaxIter(int max_iter)
Definition: hypre.cpp:4137
void SetPolyOptions(int poly_order, double poly_fraction, int eig_est_cg_iter=10)
Set parameters for polynomial smoothing.
Definition: hypre.cpp:3236
HypreDiagScale(const HypreParMatrix &A)
Definition: hypre.hpp:1294
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
Definition: hypre.cpp:3230
void BooleanMult(int alpha, const int *x, int beta, int *y)
The &quot;Boolean&quot; analog of y = alpha * A * x + beta * y, where elements in the sparsity pattern of the m...
Definition: hypre.hpp:676
void SetInterpolation(int interp_type)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1542
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:148
HypreParMatrix * EliminateCols(const Array< int > &cols)
Definition: hypre.cpp:2265
HypreParVector * X
Definition: hypre.hpp:1042
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1666
Type
Enumeration defining IDs for some classes derived from Operator.
Definition: operator.hpp:255
static MemoryType GetHostMemoryType()
Get the current Host MemoryType. This is the MemoryType used by most MFEM classes when allocating mem...
Definition: device.hpp:264
PCG solver in hypre.
Definition: hypre.hpp:1116
GMRES solver in hypre.
Definition: hypre.hpp:1175
signed char OwnsOffd() const
Get offd ownership flag.
Definition: hypre.hpp:545
void HostRead() const
Update the internal hypre_ParCSRMatrix object, A, to be on host.
Definition: hypre.hpp:782
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
signed char OwnsColMap() const
Get colmap ownership flag.
Definition: hypre.hpp:547
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:5209
Memory< double > auxX
Definition: hypre.hpp:1044
HypreParVector * X
Definition: hypre.hpp:902
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:5986
HYPRE_Int HYPRE_BigInt
void SetData(double *data_)
Sets the data of the Vector and the hypre_ParVector to data_.
Definition: hypre.cpp:309
void SetNumModes(int num_eigs)
Definition: hypre.cpp:5919
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1341
HypreParMatrix(hypre_ParCSRMatrix *a, bool owner=true)
Converts hypre&#39;s format to HypreParMatrix.
Definition: hypre.hpp:420
virtual ~HypreDiagScale()
Definition: hypre.hpp:1310
const HypreParVector & GetEigenvector(unsigned int i) const
Extract a single eigenvector.
Definition: hypre.cpp:6022
int GetNumRows() const
Returns the number of rows in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:615
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1102
const HYPRE_BigInt * ColPart() const
Returns the column partitioning (const version)
Definition: hypre.hpp:575
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.hpp:1065
HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A + B.
Definition: hypre.cpp:2613
double max_eig_est
Maximal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:940
A class to initialize the size of a Tensor.
Definition: dtensor.hpp:54
void SetRelTol(double rel_tol)
Definition: hypre.cpp:5613
HypreAME(MPI_Comm comm)
Definition: hypre.cpp:5877
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A^t * x + beta * y.
Definition: hypre.cpp:1837
void GetNumIterations(int &num_iterations) const
Definition: hypre.hpp:1557
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1339
virtual ~HypreIdentity()
Definition: hypre.hpp:1286
void SetMaxLevels(int max_levels)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1530
void SetMaxIter(int max_iter)
Definition: hypre.cpp:5943
void SetFIRCoefficients(double max_eig)
Compute window and Chebyshev coefficients for given polynomial order.
Definition: hypre.cpp:3386
void SetRandomSeed(int s)
Definition: hypre.hpp:1792
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:2279
signed char OwnsDiag() const
Get diag ownership flag.
Definition: hypre.hpp:543
double * fir_coeffs
Combined coefficients for windowing and Chebyshev polynomials.
Definition: hypre.hpp:947
void SetIsTriangular()
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1508
void Threshold(double threshold=0.0)
Remove values smaller in absolute value than some threshold.
Definition: hypre.cpp:2109
double InnerProduct(HypreParVector *x, HypreParVector *y)
Definition: hypre.cpp:426
void Print(const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0) const
Prints the locally owned rows in parallel.
Definition: hypre.cpp:2350
virtual ~HypreSmoother()
Definition: hypre.cpp:3577
void SetPrintLevel(HYPRE_Int print_level)
Set the print level: 0 = none, 1 = setup, 2 = solve, 3 = setup+solve.
Definition: hypre.cpp:4481
Host memory; using new[] and delete[].
void SetSubSpaceProjector(Operator &proj)
Definition: hypre.hpp:1799
virtual MFEM_DEPRECATED N_Vector ToNVector()
(DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_SERIAL. ...
Definition: vector.hpp:467
constexpr MemoryClass GetHypreMemoryClass()
The MemoryClass used by Hypre objects.
Definition: hypre.hpp:124
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:629
HYPRE_BigInt * ColPart()
Returns the column partitioning.
Definition: hypre.hpp:567
void HypreReadWrite()
Prepare the HypreParVector for read and write access in hypre&#39;s device memory space, HYPRE_MEMORY_DEVICE.
Definition: hypre.cpp:325
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:4653
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
Definition: hypre.cpp:2020
MFEM_DEPRECATED void SetZeroInintialIterate()
deprecated: use SetZeroInitialIterate()
Definition: hypre.hpp:1145
void SetTaubinOptions(double lambda, double mu, int iter)
Set parameters for Taubin&#39;s lambda-mu method.
Definition: hypre.cpp:3244
HypreParVector * B
Right-hand side and solution vectors.
Definition: hypre.hpp:902
HYPRE_BigInt * RowPart()
Returns the row partitioning.
Definition: hypre.hpp:563
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:4402
virtual void MultTranspose(const Vector &b, Vector &x) const
Apply transpose of the smoother to relax the linear system Ax=b.
Definition: hypre.cpp:3567
A simple singleton class for hypre&#39;s global settings, that 1) calls HYPRE_Init() and sets some GPU-re...
Definition: hypre.hpp:59
int dim
Definition: ex24.cpp:53
virtual ~HyprePCG()
Definition: hypre.cpp:3900
int GetNumCols() const
Returns the number of columns in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:622
void SetZeroInitialIterate()
non-hypre setting
Definition: hypre.hpp:1148
void SetOperator(Operator &A)
Definition: hypre.cpp:5653
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:3800
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1100
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1380
int relax_times
Number of relaxation sweeps.
Definition: hypre.hpp:916
void GetEigenvalues(Array< double > &eigenvalues) const
Collect the converged eigenvalues.
Definition: hypre.cpp:5998
void SetCycleType(int cycle_type)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1554
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
Definition: hypre.cpp:2066
void SetPrintLevel(int logging)
Definition: hypre.cpp:5949
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:437
int eig_est_cg_iter
Number of CG iterations to determine eigenvalue estimates.
Definition: hypre.hpp:938
void SetOwnership(int own)
Sets ownership of the internal hypre_ParVector.
Definition: hypre.hpp:238
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:1026
Flexible GMRES solver in hypre.
Definition: hypre.hpp:1226
ErrorMode error_mode
How to treat hypre errors.
Definition: hypre.hpp:1050
HYPRE_BigInt NNZ() const
Returns the global number of nonzeros.
Definition: hypre.hpp:559
const HypreParVector & GetEigenvector(unsigned int i) const
Extract a single eigenvector.
Definition: hypre.cpp:5722
HYPRE_BigInt * GetColStarts() const
Return the parallel column partitioning array.
Definition: hypre.hpp:644
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
ILU Setup function.
Definition: hypre.hpp:1433
const HYPRE_BigInt * Partitioning() const
Returns the parallel row/column partitioning.
Definition: hypre.hpp:217
HypreParVector & operator=(double d)
Set constant values.
Definition: hypre.cpp:278
const double alpha
Definition: ex15.cpp:369
Vector data type.
Definition: vector.hpp:60
void SetAbsTol(double tol)
Definition: hypre.cpp:3965
void PrintCommPkg(std::ostream &out=mfem::out) const
Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
Definition: hypre.cpp:2394
HypreParMatrix & operator=(double value)
Initialize all entries with value.
Definition: hypre.hpp:695
HypreParMatrix * ExtractSubmatrix(const Array< int > &indices, double threshold=0.0) const
Definition: hypre.cpp:1626
hypre_ParCSRMatrix * StealData()
Changes the ownership of the matrix.
Definition: hypre.cpp:1363
void SetPrecondUsageMode(int pcg_mode)
Definition: hypre.cpp:5638
void SetOperator(const HypreParMatrix &A)
Definition: hypre.cpp:5964
void SetSystemsOptions(int dim, bool order_bynodes=false)
Definition: hypre.cpp:4672
virtual ~HypreParMatrix()
Calls hypre&#39;s destroy function.
Definition: hypre.hpp:833
ID for class HypreParMatrix.
Definition: operator.hpp:259
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A * x + beta * y.
Definition: hypre.cpp:1703
void DropSmallEntries(double tol)
Wrapper for hypre_ParCSRMatrixDropSmallEntries in different versions of hypre. Drop off-diagonal entr...
Definition: hypre.cpp:2194
void SetAggressiveCoarsening(int num_levels)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1572
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3970
RefCoord s[3]
const HYPRE_BigInt * RowPart() const
Returns the row partitioning (const version)
Definition: hypre.hpp:571
double omega
SOR parameter (usually in (0,2))
Definition: hypre.hpp:920
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_BigInt *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:1879
Base class for solvers.
Definition: operator.hpp:651
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1626
void SetNodal(int blocksize)
Expert option - consult hypre documentation/team.
Definition: hypre.hpp:1565
int GetOwnership() const
Gets ownership of the internal hypre_ParVector.
Definition: hypre.hpp:241
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
PCG Setup function.
Definition: hypre.hpp:1161
void WrapHypreParCSRMatrix(hypre_ParCSRMatrix *a, bool owner=true)
Converts hypre&#39;s format to HypreParMatrix.
Definition: hypre.cpp:620
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:337
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:1378
void HostReadWrite()
Update the internal hypre_ParCSRMatrix object, A, to be on host.
Definition: hypre.hpp:788
void AbsMult(double a, const Vector &x, double b, Vector &y) const
Computes y = a * |A| * x + b * y, using entry-wise absolute values of the matrix A.
Definition: hypre.cpp:1845
Memory< double > auxB
Auxiliary buffers for the case when the input or output arrays in methods like Mult(const Vector &amp;...
Definition: hypre.hpp:906
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1580
~HypreParVector()
Calls hypre&#39;s destroy function.
Definition: hypre.cpp:408
void SetNumModes(int num_eigs)
Definition: hypre.hpp:1790
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:3823
HypreParMatrix * HypreParMatrixFromBlocks(Array2D< HypreParMatrix * > &blocks, Array2D< double > *blockCoeff)
Returns a merged hypre matrix constructed from hypre matrix blocks.
Definition: hypre.cpp:2822
HypreParVector * Z
Definition: hypre.hpp:908
void SetSingularProblem()
Set this option when solving a curl-curl problem with zero mass term.
Definition: hypre.hpp:1621
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:5767
virtual ~HypreFGMRES()
Definition: hypre.cpp:4238
HypreParVector * V
Temporary vectors.
Definition: hypre.hpp:908
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:3224
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:438
MemoryClass
Memory classes identify sets of memory types.
Definition: mem_manager.hpp:73
void WrapMemoryReadWrite(Memory< double > &mem)
Replace the HypreParVector&#39;s data with the given Memory, mem, and prepare the vector for read and wri...
Definition: hypre.cpp:358
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
GMRES Setup function.
Definition: hypre.hpp:1212
int poly_order
Order of the smoothing polynomial.
Definition: hypre.hpp:922
double lambda
Taubin&#39;s lambda-mu method parameters.
Definition: hypre.hpp:929
void WrapMemoryRead(const Memory< double > &mem)
Replace the HypreParVector&#39;s data with the given Memory, mem, and prepare the vector for read access ...
Definition: hypre.cpp:343
MFEM_DEPRECATED void SetZeroInintialIterate()
deprecated: use SetZeroInitialIterate()
Definition: hypre.hpp:1253
MFEM_DEPRECATED HypreParMatrix * GetData()
Deprecated. Use HypreDiagScale::GetData() const instead.
Definition: hypre.hpp:1307
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:1668
void HypreReadWrite()
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
Definition: hypre.hpp:806