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