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