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