MFEM  v3.1
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 
30 #ifdef HYPRE_COMPLEX
31 #error "MFEM does not work with HYPRE's complex numbers support"
32 #endif
33 
34 #include "sparsemat.hpp"
35 #include "hypre_parcsr.hpp"
36 
37 namespace mfem
38 {
39 
40 class ParFiniteElementSpace;
41 class HypreParMatrix;
42 
43 namespace internal
44 {
45 
46 // Convert a HYPRE_Int to int
47 inline int to_int(HYPRE_Int i)
48 {
49 #ifdef HYPRE_BIGINT
50  MFEM_ASSERT(HYPRE_Int(int(i)) == i, "overflow converting HYPRE_Int to int");
51 #endif
52  return int(i);
53 }
54 
55 }
56 
58 class HypreParVector : public Vector
59 {
60 private:
61  int own_ParVector;
62 
64  hypre_ParVector *x;
65 
66  friend class HypreParMatrix;
67 
68  // Set Vector::data and Vector::size from *x
69  inline void _SetDataAndSize_();
70 
71 public:
74  HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col);
77  HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, double *_data,
78  HYPRE_Int *col);
82  HypreParVector(HypreParMatrix &A, int tr = 0);
84  HypreParVector(HYPRE_ParVector y);
87 
89  MPI_Comm GetComm() { return x->comm; }
90 
92  inline HYPRE_Int *Partitioning() { return x->partitioning; }
93 
95  inline HYPRE_Int GlobalSize() { return x->global_size; }
96 
98  operator hypre_ParVector*() const;
99 #ifndef HYPRE_PAR_VECTOR_STRUCT
100  operator HYPRE_ParVector() const;
102 #endif
103  hypre_ParVector *StealParVector() { own_ParVector = 0; return x; }
105 
107  void SetOwnership(int own) { own_ParVector = own; }
108 
110  int GetOwnership() const { return own_ParVector; }
111 
113  Vector* GlobalVector() const;
114 
116  HypreParVector& operator= (double d);
119 
125  void SetData(double *_data);
126 
128  HYPRE_Int Randomize(HYPRE_Int seed);
129 
131  void Print(const char *fname) const;
132 
134  ~HypreParVector();
135 };
136 
138 double InnerProduct(HypreParVector &x, HypreParVector &y);
139 double InnerProduct(HypreParVector *x, HypreParVector *y);
140 
141 
143 class HypreParMatrix : public Operator
144 {
145 private:
147  hypre_ParCSRMatrix *A;
148 
150  mutable HypreParVector *X, *Y;
151 
152  // Flags indicating ownership of A->diag->{i,j,data}, A->offd->{i,j,data},
153  // and A->col_map_offd.
154  // The possible values for diagOwner are:
155  // -1: no special treatment of A->diag (default)
156  // 0: prevent hypre from destroying A->diag->{i,j,data}
157  // 1: same as 0, plus take ownership of A->diag->{i,j}
158  // 2: same as 0, plus take ownership of A->diag->data
159  // 3: same as 0, plus take ownership of A->diag->{i,j,data}
160  // The same values and rules apply to offdOwner and A->offd.
161  // The possible values for colMapOwner are:
162  // -1: no special treatment of A->col_map_offd (default)
163  // 0: prevent hypre from destroying A->col_map_offd
164  // 1: same as 0, plus take ownership of A->col_map_offd
165  // All owned arrays are destroyed with 'delete []'.
166  char diagOwner, offdOwner, colMapOwner;
167 
168  // Does the object own the pointer A?
169  char ParCSROwner;
170 
171  // Initialize with defaults. Does not initialize inherited members.
172  void Init();
173 
174  // Delete all owned data. Does not perform re-initialization with defaults.
175  void Destroy();
176 
177  // Copy (shallow or deep, based on HYPRE_BIGINT) the I and J arrays from csr
178  // to hypre_csr. Shallow copy the data. Return the appropriate ownership
179  // flag.
180  static char CopyCSR(SparseMatrix *csr, hypre_CSRMatrix *hypre_csr);
181  // Copy (shallow or deep, based on HYPRE_BIGINT) the I and J arrays from
182  // bool_csr to hypre_csr. Allocate the data array and set it to all ones.
183  // Return the appropriate ownership flag.
184  static char CopyBoolCSR(Table *bool_csr, hypre_CSRMatrix *hypre_csr);
185 
186  // Copy the j array of a hypre_CSRMatrix to the given J array, converting
187  // the indices from HYPRE_Int to int.
188  static void CopyCSR_J(hypre_CSRMatrix *hypre_csr, int *J);
189 
190 public:
192  HypreParMatrix();
193 
195  HypreParMatrix(hypre_ParCSRMatrix *a)
196  {
197  Init();
198  A = a;
199  height = GetNumRows();
200  width = GetNumCols();
201  }
202 
206  HypreParMatrix(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *row_starts,
207  SparseMatrix *diag);
208 
212  HypreParMatrix(MPI_Comm comm, HYPRE_Int global_num_rows,
213  HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
214  HYPRE_Int *col_starts, SparseMatrix *diag);
215 
218  HypreParMatrix(MPI_Comm comm, HYPRE_Int global_num_rows,
219  HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
220  HYPRE_Int *col_starts, SparseMatrix *diag, SparseMatrix *offd,
221  HYPRE_Int *cmap);
222 
225  HypreParMatrix(MPI_Comm comm,
226  HYPRE_Int global_num_rows, HYPRE_Int global_num_cols,
227  HYPRE_Int *row_starts, HYPRE_Int *col_starts,
228  HYPRE_Int *diag_i, HYPRE_Int *diag_j, double *diag_data,
229  HYPRE_Int *offd_i, HYPRE_Int *offd_j, double *offd_data,
230  HYPRE_Int offd_num_cols, HYPRE_Int *offd_col_map);
231 
233  HypreParMatrix(MPI_Comm comm, HYPRE_Int *row_starts, HYPRE_Int *col_starts,
234  SparseMatrix *a);
235 
238  HypreParMatrix(MPI_Comm comm, HYPRE_Int global_num_rows,
239  HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
240  HYPRE_Int *col_starts, Table *diag);
241 
245  HypreParMatrix(MPI_Comm comm, int id, int np, HYPRE_Int *row, HYPRE_Int *col,
246  HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd,
247  HYPRE_Int *j_offd, HYPRE_Int *cmap, HYPRE_Int cmap_size);
248 
253  HypreParMatrix(MPI_Comm comm, int nrows, HYPRE_Int glob_nrows,
254  HYPRE_Int glob_ncols, int *I, HYPRE_Int *J,
255  double *data, HYPRE_Int *rows, HYPRE_Int *cols);
256 
258  void MakeRef(const HypreParMatrix &master);
259 
261  MPI_Comm GetComm() const { return A->comm; }
262 
264  operator hypre_ParCSRMatrix*() { return A; }
265 #ifndef HYPRE_PAR_CSR_MATRIX_STRUCT
266  operator HYPRE_ParCSRMatrix() { return (HYPRE_ParCSRMatrix) A; }
268 #endif
269  hypre_ParCSRMatrix* StealData();
271 
273  void SetOwnerFlags(char diag, char offd, char colmap)
274  { diagOwner = diag, offdOwner = offd, colMapOwner = colmap; }
275 
277  char OwnsDiag() const { return diagOwner; }
279  char OwnsOffd() const { return offdOwner; }
281  char OwnsColMap() const { return colMapOwner; }
282 
286  void CopyRowStarts();
290  void CopyColStarts();
291 
293  inline HYPRE_Int NNZ() { return A->num_nonzeros; }
295  inline HYPRE_Int *RowPart() { return A->row_starts; }
297  inline HYPRE_Int *ColPart() { return A->col_starts; }
299  inline HYPRE_Int M() { return A->global_num_rows; }
301  inline HYPRE_Int N() { return A->global_num_cols; }
302 
304  void GetDiag(Vector &diag) const;
306  void GetDiag(SparseMatrix &diag) const;
308  void GetOffd(SparseMatrix &offd, HYPRE_Int* &cmap) const;
309 
312  void GetBlocks(Array2D<HypreParMatrix*> &blocks,
313  bool interleaved_rows = false,
314  bool interleaved_cols = false) const;
315 
318 
320  int GetNumRows() const
321  {
322  return internal::to_int(
323  hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(A)));
324  }
325 
327  int GetNumCols() const
328  {
329  return internal::to_int(
330  hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(A)));
331  }
332 
333  HYPRE_Int GetGlobalNumRows() const
334  { return hypre_ParCSRMatrixGlobalNumRows(A); }
335 
336  HYPRE_Int GetGlobalNumCols() const
337  { return hypre_ParCSRMatrixGlobalNumCols(A); }
338 
339  HYPRE_Int *GetRowStarts() const { return hypre_ParCSRMatrixRowStarts(A); }
340 
341  HYPRE_Int *GetColStarts() const { return hypre_ParCSRMatrixColStarts(A); }
342 
344  HYPRE_Int Mult(HypreParVector &x, HypreParVector &y,
345  double alpha = 1.0, double beta = 0.0);
347  HYPRE_Int Mult(HYPRE_ParVector x, HYPRE_ParVector y,
348  double alpha = 1.0, double beta = 0.0);
351  double alpha = 1.0, double beta = 0.0);
352 
353  void Mult(double a, const Vector &x, double b, Vector &y) const;
354  void MultTranspose(double a, const Vector &x, double b, Vector &y) const;
355 
356  virtual void Mult(const Vector &x, Vector &y) const
357  { Mult(1.0, x, 0.0, y); }
358  virtual void MultTranspose(const Vector &x, Vector &y) const
359  { MultTranspose(1.0, x, 0.0, y); }
360 
363  void BooleanMult(int alpha, int *x, int beta, int *y)
364  {
365  internal::hypre_ParCSRMatrixBooleanMatvec(A, alpha, x, beta, y);
366  }
367 
376  HYPRE_Int* row_starts = NULL) const;
377 
379  void ScaleRows(const Vector & s);
381  void InvScaleRows(const Vector & s);
383  void operator*=(double s);
384 
386  void Threshold(double threshold = 0.0);
387 
389  void EliminateZeroRows() { hypre_ParCSRMatrixFixZeroRows(A); }
390 
393  void EliminateRowsCols(const Array<int> &rows_cols, const HypreParVector &X,
394  HypreParVector &B);
395 
399  HypreParMatrix* EliminateRowsCols(const Array<int> &rows_cols);
400 
402  void Print(const char *fname, HYPRE_Int offi = 0, HYPRE_Int offj = 0);
404  void Read(MPI_Comm comm, const char *fname);
405 
407  virtual ~HypreParMatrix() { Destroy(); }
408 };
409 
411 HypreParMatrix * ParMult(HypreParMatrix *A, HypreParMatrix *B);
412 
414 HypreParMatrix * RAP(HypreParMatrix *A, HypreParMatrix *P);
416 HypreParMatrix * RAP(HypreParMatrix * Rt, HypreParMatrix *A, HypreParMatrix *P);
417 
421 void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae,
422  const Array<int> &ess_dof_list, const Vector &X, Vector &B);
423 
424 
426 class HypreSmoother : public Solver
427 {
428 protected:
432  mutable HypreParVector *B, *X;
434  mutable HypreParVector *V, *Z;
436  mutable HypreParVector *X0, *X1;
437 
440  int type;
444  double relax_weight;
446  double omega;
453 
455  double lambda;
456  double mu;
458 
460  double *l1_norms;
462  double max_eig_est;
464  double min_eig_est;
466  double window_params[3];
467 
469  double* fir_coeffs;
470 
471 public:
482  enum Type { Jacobi = 0, l1Jacobi = 1, l1GS = 2, l1GStr = 4, lumpedJacobi = 5,
483  GS = 6, Chebyshev = 16, Taubin = 1001, FIR = 1002
484  };
485 
486  HypreSmoother();
487 
489  int relax_times = 1, double relax_weight = 1.0,
490  double omega = 1.0, int poly_order = 2,
491  double poly_fraction = .3);
492 
496  void SetSOROptions(double relax_weight, double omega);
498  void SetPolyOptions(int poly_order, double poly_fraction);
500  void SetTaubinOptions(double lambda, double mu, int iter);
501 
503  void SetWindowByName(const char* window_name);
505  void SetWindowParameters(double a, double b, double c);
507  void SetFIRCoefficients(double max_eig);
508 
511  virtual void SetOperator(const Operator &op);
512 
514  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
515  virtual void Mult(const Vector &b, Vector &x) const;
516 
517  virtual ~HypreSmoother();
518 };
519 
520 
522 class HypreSolver : public Solver
523 {
524 protected:
527 
529  mutable HypreParVector *B, *X;
530 
532  mutable int setup_called;
533 
534 public:
535  HypreSolver();
536 
538 
540  virtual operator HYPRE_Solver() const = 0;
541 
543  virtual HYPRE_PtrToParSolverFcn SetupFcn() const = 0;
545  virtual HYPRE_PtrToParSolverFcn SolveFcn() const = 0;
546 
547  virtual void SetOperator(const Operator &op)
548  { mfem_error("HypreSolvers do not support SetOperator!"); }
549 
551  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
552  virtual void Mult(const Vector &b, Vector &x) const;
553 
554  virtual ~HypreSolver();
555 };
556 
558 class HyprePCG : public HypreSolver
559 {
560 private:
561  HYPRE_Solver pcg_solver;
562 
563 public:
565 
566  void SetTol(double tol);
567  void SetMaxIter(int max_iter);
568  void SetLogging(int logging);
569  void SetPrintLevel(int print_lvl);
570 
572  void SetPreconditioner(HypreSolver &precond);
573 
577  void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0);
578 
581 
582  void GetNumIterations(int &num_iterations)
583  {
584  HYPRE_Int num_it;
585  HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_it);
586  num_iterations = internal::to_int(num_it);
587  }
588 
590  virtual operator HYPRE_Solver() const { return pcg_solver; }
591 
593  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
594  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSetup; }
596  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
597  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSolve; }
598 
600  virtual void Mult(const HypreParVector &b, HypreParVector &x) const;
601  using HypreSolver::Mult;
602 
603  virtual ~HyprePCG();
604 };
605 
607 class HypreGMRES : public HypreSolver
608 {
609 private:
610  HYPRE_Solver gmres_solver;
611 
612 public:
614 
615  void SetTol(double tol);
616  void SetMaxIter(int max_iter);
617  void SetKDim(int dim);
618  void SetLogging(int logging);
619  void SetPrintLevel(int print_lvl);
620 
622  void SetPreconditioner(HypreSolver &precond);
623 
626 
628  virtual operator HYPRE_Solver() const { return gmres_solver; }
629 
631  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
632  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSetup; }
634  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
635  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSolve; }
636 
638  virtual void Mult (const HypreParVector &b, HypreParVector &x) const;
639  using HypreSolver::Mult;
640 
641  virtual ~HypreGMRES();
642 };
643 
646 {
647 public:
648  virtual operator HYPRE_Solver() const { return NULL; }
649 
650  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
651  { return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentitySetup; }
652  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
653  { return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentity; }
654 
655  virtual ~HypreIdentity() { }
656 };
657 
660 {
661 public:
664  virtual operator HYPRE_Solver() const { return NULL; }
665 
666  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
667  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScaleSetup; }
668  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
669  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScale; }
670 
671  HypreParMatrix* GetData() { return A; }
672  virtual ~HypreDiagScale() { }
673 };
674 
677 {
678 private:
679  HYPRE_Solver sai_precond;
680 
681 public:
683 
684  void SetSymmetry(int sym);
685 
687  virtual operator HYPRE_Solver() const { return sai_precond; }
688 
689  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
690  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSetup; }
691  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
692  { return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSolve; }
693 
694  virtual ~HypreParaSails();
695 };
696 
699 {
700 private:
701  HYPRE_Solver amg_precond;
702 
705 
707  ParFiniteElementSpace *fespace;
708 
710  void RecomputeRBMs();
711 
713  void SetDefaultOptions();
714 
715  // If amg_precond is NULL, allocates it and sets default options.
716  // Otherwise saves the options from amg_precond, destroys it, allocates a new
717  // one, and sets its options to the saved values.
718  void ResetAMGPrecond();
719 
720 public:
721  HypreBoomerAMG();
722 
724 
725  virtual void SetOperator(const Operator &op);
726 
730  void SetSystemsOptions(int dim);
731 
738 
739  void SetPrintLevel(int print_level)
740  { HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level); }
741 
743  virtual operator HYPRE_Solver() const { return amg_precond; }
744 
745  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
746  { return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSetup; }
747  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
748  { return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSolve; }
749 
750  virtual ~HypreBoomerAMG();
751 };
752 
754 HypreParMatrix* DiscreteGrad(ParFiniteElementSpace *edge_fespace,
755  ParFiniteElementSpace *vert_fespace);
757 HypreParMatrix* DiscreteCurl(ParFiniteElementSpace *face_fespace,
758  ParFiniteElementSpace *edge_fespace);
759 
761 class HypreAMS : public HypreSolver
762 {
763 private:
764  HYPRE_Solver ams;
765 
767  HypreParVector *x, *y, *z;
769  HypreParMatrix *G;
771  HypreParMatrix *Pi, *Pix, *Piy, *Piz;
772 
773 public:
775 
776  void SetPrintLevel(int print_lvl);
777 
779  void SetSingularProblem() { HYPRE_AMSSetBetaPoissonMatrix(ams, NULL); }
780 
782  virtual operator HYPRE_Solver() const { return ams; }
783 
784  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
785  { return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSetup; }
786  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
787  { return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSolve; }
788 
789  virtual ~HypreAMS();
790 };
791 
793 class HypreADS : public HypreSolver
794 {
795 private:
796  HYPRE_Solver ads;
797 
799  HypreParVector *x, *y, *z;
801  HypreParMatrix *G;
803  HypreParMatrix *C;
805  HypreParMatrix *ND_Pi, *ND_Pix, *ND_Piy, *ND_Piz;
807  HypreParMatrix *RT_Pi, *RT_Pix, *RT_Piy, *RT_Piz;
808 
809 public:
811 
812  void SetPrintLevel(int print_lvl);
813 
815  virtual operator HYPRE_Solver() const { return ads; }
816 
817  virtual HYPRE_PtrToParSolverFcn SetupFcn() const
818  { return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSetup; }
819  virtual HYPRE_PtrToParSolverFcn SolveFcn() const
820  { return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSolve; }
821 
822  virtual ~HypreADS();
823 };
824 
848 {
849 private:
850  MPI_Comm comm;
851  int myid;
852  int numProcs;
853  int nev; // Number of desired eigenmodes
854  int seed; // Random seed used for initial vectors
855 
856  HYPRE_Int glbSize; // Global number of DoFs in the linear system
857  HYPRE_Int * part; // Row partitioning of the linear system
858 
859  // Pointer to HYPRE's solver struct
860  HYPRE_Solver lobpcg_solver;
861 
862  // Interface for matrix storage type
863  mv_InterfaceInterpreter interpreter;
864 
865  // Interface for setting up and performing matrix-vector products
866  HYPRE_MatvecFunctions matvec_fn;
867 
868  // Eigenvalues
869  Array<double> eigenvalues;
870 
871  // Forward declaration
872  class HypreMultiVector;
873 
874  // MultiVector to store eigenvectors
875  HypreMultiVector * multi_vec;
876 
877  // Empty vectors used to setup the matrices and preconditioner
878  HypreParVector * x;
879 
881  class HypreMultiVector
882  {
883  private:
884  // Pointer to hypre's multi-vector object
885  mv_MultiVectorPtr mv_ptr;
886 
887  // Wrappers for each member of the multivector
888  HypreParVector ** hpv;
889 
890  // Number of vectors in the multivector
891  int nv;
892 
893  public:
894  HypreMultiVector(int n, HypreParVector & v,
895  mv_InterfaceInterpreter & interpreter);
896  ~HypreMultiVector();
897 
899  void Randomize(HYPRE_Int seed);
900 
902  HypreParVector & GetVector(unsigned int i);
903 
905  HypreParVector ** StealVectors();
906 
907  operator mv_MultiVectorPtr() const { return mv_ptr; }
908 
909  mv_MultiVectorPtr & GetMultiVector() { return mv_ptr; }
910  };
911 
912  static void * OperatorMatvecCreate( void *A, void *x );
913  static HYPRE_Int OperatorMatvec( void *matvec_data,
914  HYPRE_Complex alpha,
915  void *A,
916  void *x,
917  HYPRE_Complex beta,
918  void *y );
919  static HYPRE_Int OperatorMatvecDestroy( void *matvec_data );
920 
921  static HYPRE_Int PrecondSolve(void *solver,
922  void *A,
923  void *b,
924  void *x);
925  static HYPRE_Int PrecondSetup(void *solver,
926  void *A,
927  void *b,
928  void *x);
929 
930 public:
931  HypreLOBPCG(MPI_Comm comm);
932  ~HypreLOBPCG();
933 
934  void SetTol(double tol);
935  void SetMaxIter(int max_iter);
936  void SetPrintLevel(int logging);
937  void SetNumModes(int num_eigs) { nev = num_eigs; }
938  void SetPrecondUsageMode(int pcg_mode);
939  void SetRandomSeed(int s) { seed = s; }
940 
941  // The following four methods support general operators
942  void SetPreconditioner(Solver & precond);
943  void SetOperator(Operator & A);
944  void SetMassMatrix(Operator & M);
945 
947  void Solve();
948 
950  void GetEigenvalues(Array<double> & eigenvalues);
951 
953  HypreParVector & GetEigenvector(unsigned int i);
954 
956  HypreParVector ** StealEigenvectors() { return multi_vec->StealVectors(); }
957 };
958 
981 class HypreAME
982 {
983 private:
984  int myid;
985  int numProcs;
986  int nev; // Number of desired eigenmodes
987  bool setT;
988 
989  // Pointer to HYPRE's AME solver struct
990  HYPRE_Solver ame_solver;
991 
992  // Pointer to HYPRE's AMS solver struct
993  HypreSolver * ams_precond;
994 
995  // Eigenvalues
996  HYPRE_Real * eigenvalues;
997 
998  // MultiVector to store eigenvectors
999  HYPRE_ParVector * multi_vec;
1000 
1001  HypreParVector ** eigenvectors;
1002 
1003  void createDummyVectors();
1004 
1005 public:
1006  HypreAME(MPI_Comm comm);
1007  ~HypreAME();
1008 
1009  void SetTol(double tol);
1010  void SetMaxIter(int max_iter);
1011  void SetPrintLevel(int logging);
1012  void SetNumModes(int num_eigs);
1013 
1014  // The following four methods support operators of type HypreParMatrix.
1015  void SetPreconditioner(HypreSolver & precond);
1016  void SetOperator(HypreParMatrix & A);
1017  void SetMassMatrix(HypreParMatrix & M);
1018 
1020  void Solve();
1021 
1023  void GetEigenvalues(Array<double> & eigenvalues);
1024 
1026  HypreParVector & GetEigenvector(unsigned int i);
1027 
1030 };
1031 
1032 }
1033 
1034 #endif // MFEM_USE_MPI
1035 
1036 #endif
virtual ~HypreBoomerAMG()
Definition: hypre.cpp:2444
int to_int(string str)
void SetTol(double tol)
Definition: hypre.cpp:1917
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2075
HYPRE_Int NNZ()
Returns the global number of nonzeros.
Definition: hypre.hpp:293
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:1213
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:261
Vector * GlobalVector() const
Returns the global vector in each processor.
Definition: hypre.cpp:136
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:761
HypreParVector * X0
FIR Filter Temporary Vectors.
Definition: hypre.hpp:436
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:3102
double min_eig_est
Minimal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:464
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:652
void SetOwnerFlags(char diag, char offd, char colmap)
Explicitly set the three ownership flags, see docs for diagOwner etc.
Definition: hypre.hpp:273
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:793
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2894
void Print(const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0)
Prints the locally owned rows in parallel.
Definition: hypre.cpp:1236
MPI_Comm GetComm()
MPI communicator.
Definition: hypre.hpp:89
HYPRE_Int M()
Returns the global number of rows.
Definition: hypre.hpp:299
int setup_called
Was hypre&#39;s Setup function called already?
Definition: hypre.hpp:532
HYPRE_Int * ColPart()
Returns the column partitioning.
Definition: hypre.hpp:297
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to &#39;master&#39;.
Definition: hypre.cpp:726
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_Int *row_starts=NULL) const
Definition: hypre.cpp:969
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:963
HypreParVector * B
Right-hand side and solution vector.
Definition: hypre.hpp:529
HypreDiagScale(HypreParMatrix &A)
Definition: hypre.hpp:663
double window_params[3]
Parameters for windowing function of FIR filter.
Definition: hypre.hpp:466
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
Definition: hypre.hpp:356
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:786
HypreADS(HypreParMatrix &A, ParFiniteElementSpace *face_fespace)
Definition: hypre.cpp:2654
void SetWindowByName(const char *window_name)
Convenience function for setting canonical windowing parameters.
Definition: hypre.cpp:1596
HypreParMatrix(hypre_ParCSRMatrix *a)
Converts hypre&#39;s format to HypreParMatrix.
Definition: hypre.hpp:195
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:3307
void SetPreconditioner(HypreSolver &precond)
Definition: hypre.cpp:3273
void SetTol(double tol)
Definition: hypre.cpp:3252
virtual ~HypreGMRES()
Definition: hypre.cpp:2151
void SetPreconditioner(Solver &precond)
Definition: hypre.cpp:3025
void SetMassMatrix(Operator &M)
Definition: hypre.cpp:3080
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetPrintLevel(int logging)
Definition: hypre.cpp:3010
bool iterative_mode
If true, use the second argument of Mult as an initial guess.
Definition: operator.hpp:106
virtual ~HypreAMS()
Definition: hypre.cpp:2634
void SetTol(double tol)
Definition: hypre.cpp:2050
void SetOperator(HypreParMatrix &A)
Definition: hypre.cpp:3279
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
GMRES Solve function.
Definition: hypre.hpp:634
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.hpp:956
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:903
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2070
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
Definition: hypre.cpp:1046
double poly_fraction
Fraction of spectrum to smooth for polynomial relaxation.
Definition: hypre.hpp:450
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:1932
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:750
int poly_scale
Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}.
Definition: hypre.hpp:452
void SetTol(double tol)
Definition: hypre.cpp:2998
HypreParVector * X1
Definition: hypre.hpp:436
The identity operator as a hypre solver.
Definition: hypre.hpp:645
void SetWindowParameters(double a, double b, double c)
Set parameters for windowing function for FIR smoother.
Definition: hypre.cpp:1611
void EliminateZeroRows()
If a row contains only zeros, set its diagonal to 1.
Definition: hypre.hpp:389
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s GMRES.
Definition: hypre.cpp:2083
void Print(const char *fname) const
Prints the locally owned rows in parallel.
Definition: hypre.cpp:178
void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
Definition: hypre.cpp:1368
void SetSystemsOptions(int dim)
Definition: hypre.cpp:2326
HypreLOBPCG(MPI_Comm comm)
Definition: hypre.cpp:2969
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
Definition: hypre.cpp:1728
hypre_ParVector * StealParVector()
Changes the ownership of the the vector.
Definition: hypre.hpp:104
The BoomerAMG solver in hypre.
Definition: hypre.hpp:698
HYPRE_Int GetGlobalNumRows() const
Definition: hypre.hpp:333
void SetSymmetry(int sym)
Definition: hypre.cpp:2180
int dim
Definition: ex3.cpp:48
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator.
Definition: hypre.hpp:358
virtual ~HypreParaSails()
Definition: hypre.cpp:2185
The ParaSails preconditioner in hypre.
Definition: hypre.hpp:676
HYPRE_Int GetGlobalNumCols() const
Definition: hypre.hpp:336
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3004
double * l1_norms
l1 norms of the rows of A
Definition: hypre.hpp:460
Data type sparse matrix.
Definition: sparsemat.hpp:38
void SetLogging(int logging)
Definition: hypre.cpp:1927
void SetZeroInintialIterate()
non-hypre setting
Definition: hypre.hpp:580
Jacobi preconditioner in hypre.
Definition: hypre.hpp:659
virtual ~HypreSolver()
Definition: hypre.cpp:1899
HyprePCG(HypreParMatrix &_A)
Definition: hypre.cpp:1906
virtual HYPRE_PtrToParSolverFcn SolveFcn() const =0
hypre&#39;s internal Solve function
void SetKDim(int dim)
Definition: hypre.cpp:2060
void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0)
Definition: hypre.cpp:1945
void SetPrintLevel(int print_level)
Definition: hypre.hpp:739
void GetBlocks(Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
Definition: hypre.cpp:870
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
Definition: hypre.cpp:1852
HypreParMatrix * RAP(HypreParMatrix *A, HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:1318
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
PCG Solve function.
Definition: hypre.hpp:596
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
Definition: hypre.cpp:3354
HypreParMatrix * DiscreteCurl(ParFiniteElementSpace *face_fespace, ParFiniteElementSpace *edge_fespace)
Compute the discrete curl matrix between the ND1 and RT0 spaces.
HYPRE_Int Randomize(HYPRE_Int seed)
Set random values.
Definition: hypre.cpp:173
double relax_weight
Damping coefficient (usually &lt;= 1)
Definition: hypre.hpp:444
Parallel smoothers in hypre.
Definition: hypre.hpp:426
void SetElasticityOptions(ParFiniteElementSpace *fespace)
Definition: hypre.cpp:2408
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:828
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:745
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:650
void BooleanMult(int alpha, int *x, int beta, int *y)
Definition: hypre.hpp:363
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:430
HypreAMS(HypreParMatrix &A, ParFiniteElementSpace *edge_fespace)
Definition: hypre.cpp:2455
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:668
virtual void SetOperator(const Operator &op)
Definition: hypre.cpp:1618
virtual HYPRE_PtrToParSolverFcn SetupFcn() const =0
hypre&#39;s internal Setup function
void SetLogging(int logging)
Definition: hypre.cpp:2065
char OwnsColMap() const
Get colmap ownership flag.
Definition: hypre.hpp:281
virtual ~HypreADS()
Definition: hypre.cpp:2872
HypreParMatrix()
An empty matrix to be used as a reference to an existing matrix.
Definition: hypre.cpp:211
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:666
void SetMaxIter(int max_iter)
Definition: hypre.cpp:1922
HypreParMatrix * GetData()
Definition: hypre.hpp:671
HYPRE_Int * RowPart()
Returns the row partitioning.
Definition: hypre.hpp:295
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
Definition: hypre.cpp:1576
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:58
HypreParVector * X
Definition: hypre.hpp:529
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:817
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Definition: hypre.cpp:3090
PCG solver in hypre.
Definition: hypre.hpp:558
GMRES solver in hypre.
Definition: hypre.hpp:607
char OwnsOffd() const
Get offd ownership flag.
Definition: hypre.hpp:279
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2649
HypreParVector * X
Definition: hypre.hpp:432
void mfem_error(const char *msg)
Definition: error.cpp:23
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:3301
void SetNumModes(int num_eigs)
Definition: hypre.cpp:3244
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:691
virtual ~HypreDiagScale()
Definition: hypre.hpp:672
int GetNumRows() const
Returns the number of rows in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:320
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.hpp:547
double max_eig_est
Maximal eigenvalue estimate for polynomial smoothing.
Definition: hypre.hpp:462
HypreAME(MPI_Comm comm)
Definition: hypre.cpp:3202
HypreParMatrix * Transpose()
Returns the transpose of *this.
Definition: hypre.cpp:892
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:689
virtual ~HypreIdentity()
Definition: hypre.hpp:655
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3258
void SetFIRCoefficients(double max_eig)
Compute window and Chebyshev coefficients for given polynomial order.
Definition: hypre.cpp:1688
void SetRandomSeed(int s)
Definition: hypre.hpp:939
double * fir_coeffs
Combined coefficients for windowing and Chebyshev polynomials.
Definition: hypre.hpp:469
void Threshold(double threshold=0.0)
Remove values smaller in absolute value than some threshold.
Definition: hypre.cpp:1160
double InnerProduct(HypreParVector *x, HypreParVector *y)
Definition: hypre.cpp:192
HYPRE_Int GlobalSize()
Returns the global number of rows.
Definition: hypre.hpp:95
virtual ~HypreSmoother()
Definition: hypre.cpp:1818
HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col)
Definition: hypre.cpp:50
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: hypre.cpp:2309
void SetData(double *_data)
Definition: hypre.cpp:168
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
Definition: hypre.cpp:1081
HYPRE_Int N()
Returns the global number of columns.
Definition: hypre.hpp:301
void SetTaubinOptions(double lambda, double mu, int iter)
Set parameters for Taubin&#39;s lambda-mu method.
Definition: hypre.cpp:1588
HypreParVector * B
Right-hand side and solution vectors.
Definition: hypre.hpp:432
void SetMassMatrix(HypreParMatrix &M)
Definition: hypre.cpp:3294
void SetZeroInintialIterate()
non-hypre setting
Definition: hypre.hpp:625
virtual ~HyprePCG()
Definition: hypre.cpp:2026
int GetNumCols() const
Returns the number of columns in the diagonal block of the ParCSRMatrix.
Definition: hypre.hpp:327
void SetOperator(Operator &A)
Definition: hypre.cpp:3034
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:1937
void SetPolyOptions(int poly_order, double poly_fraction)
Set parameters for polynomial smoothing.
Definition: hypre.cpp:1582
void GetNumIterations(int &num_iterations)
Definition: hypre.hpp:582
void CopyColStarts()
Definition: hypre.cpp:787
int relax_times
Number of relaxation sweeps.
Definition: hypre.hpp:442
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
Definition: hypre.cpp:1122
void SetPrintLevel(int logging)
Definition: hypre.cpp:3264
void SetOwnership(int own)
Sets ownership of the internal hypre_ParVector.
Definition: hypre.hpp:107
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:522
char OwnsDiag() const
Get diag ownership flag.
Definition: hypre.hpp:277
HypreParVector & operator=(double d)
Set constant values.
Definition: hypre.cpp:146
Vector data type.
Definition: vector.hpp:33
HypreGMRES(HypreParMatrix &_A)
Definition: hypre.cpp:2032
hypre_ParCSRMatrix * StealData()
Changes the ownership of the the matrix.
Definition: hypre.cpp:736
void SetPrecondUsageMode(int pcg_mode)
Definition: hypre.cpp:3019
virtual ~HypreParMatrix()
Calls hypre&#39;s destroy function.
Definition: hypre.hpp:407
HypreParMatrix * A
The linear system matrix.
Definition: hypre.hpp:526
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2055
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:864
HYPRE_Int * GetColStarts() const
Definition: hypre.hpp:341
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
Definition: hypre.cpp:3343
double omega
SOR parameter (usually in (0,2))
Definition: hypre.hpp:446
HYPRE_Int * GetRowStarts() const
Definition: hypre.hpp:339
Base class for solvers.
Definition: operator.hpp:102
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre&#39;s internal Setup function
Definition: hypre.hpp:784
int GetOwnership() const
Gets ownership of the internal hypre_ParVector.
Definition: hypre.hpp:110
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
PCG Setup function.
Definition: hypre.hpp:593
Abstract operator.
Definition: operator.hpp:21
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:143
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:747
~HypreParVector()
Calls hypre&#39;s destroy function.
Definition: hypre.cpp:183
void SetNumModes(int num_eigs)
Definition: hypre.hpp:937
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:1958
HypreParVector * Z
Definition: hypre.hpp:434
void SetSingularProblem()
Set this option when solving a curl-curl problem with zero mass term.
Definition: hypre.hpp:779
void Solve()
Solve the eigenproblem.
Definition: hypre.cpp:3108
void Read(MPI_Comm comm, const char *fname)
Reads the matrix from a file.
Definition: hypre.cpp:1241
HypreParVector * V
Temporary vectors.
Definition: hypre.hpp:434
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:1570
HypreParaSails(HypreParMatrix &A)
Definition: hypre.cpp:2157
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
GMRES Setup function.
Definition: hypre.hpp:631
int poly_order
Order of the smoothing polynomial.
Definition: hypre.hpp:448
HYPRE_Int * Partitioning()
Returns the row partitioning.
Definition: hypre.hpp:92
double lambda
Taubin&#39;s lambda-mu method parameters.
Definition: hypre.hpp:455
HypreParMatrix * ParMult(HypreParMatrix *A, HypreParMatrix *B)
Returns the matrix A * B.
Definition: hypre.cpp:1308
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre&#39;s internal Solve function
Definition: hypre.hpp:819