15 #include "../config/config.hpp" 26 #include "_hypre_parcsr_mv.h" 27 #include "_hypre_parcsr_ls.h" 28 #include "temp_multivector.h" 29 #include "../general/globals.hpp" 32 #error "MFEM does not work with HYPRE's complex numbers support" 35 #if defined(HYPRE_USING_GPU) && \ 36 !(defined(HYPRE_USING_CUDA) || defined(HYPRE_USING_HIP)) 37 #error "Unsupported GPU build of HYPRE! Only CUDA and HIP builds are supported." 39 #if defined(HYPRE_USING_CUDA) && !defined(MFEM_USE_CUDA) 40 #error "MFEM_USE_CUDA=YES is required when HYPRE is built with CUDA!" 42 #if defined(HYPRE_USING_HIP) && !defined(MFEM_USE_HIP) 43 #error "MFEM_USE_HIP=YES is required when HYPRE is built with HIP!" 49 #if defined(HYPRE_USING_CUDA) 50 #define MFEM_HYPRE_FORALL(i, N,...) CuWrap1D(N, [=] MFEM_DEVICE \ 51 (int i) {__VA_ARGS__}) 52 #elif defined(HYPRE_USING_HIP) 53 #define MFEM_HYPRE_FORALL(i, N,...) HipWrap1D(N, [=] MFEM_DEVICE \ 54 (int i) {__VA_ARGS__}) 56 #define MFEM_HYPRE_FORALL(i, N,...) for (int i = 0; i < N; i++) { __VA_ARGS__ } 65 class ParFiniteElementSpace;
80 static void Init() { Instance(); }
97 void SetDefaultOptions();
100 static Hypre &Instance()
106 bool finalized =
false;
113 template <
typename int_type>
114 inline int to_int(int_type i)
116 MFEM_ASSERT(int_type(
int(i)) == i,
"overflow converting int_type to int");
121 template <>
inline int to_int(
int i) {
return i; }
126 inline int to_int(HYPRE_Int i)
128 MFEM_ASSERT(HYPRE_Int(
int(i)) == i,
"overflow converting HYPRE_Int to int");
139 #if !defined(HYPRE_USING_GPU) 141 #elif defined(HYPRE_USING_UNIFIED_MEMORY) 151 #if !defined(HYPRE_USING_GPU) 153 #elif defined(HYPRE_USING_UNIFIED_MEMORY) 172 inline void _SetDataAndSize_();
179 own_ParVector =
false;
242 operator hypre_ParVector*()
const {
return x; }
243 #ifndef HYPRE_PAR_VECTOR_STRUCT 245 operator HYPRE_ParVector()
const {
return (HYPRE_ParVector) x; }
323 void Print(
const char *fname)
const;
326 void Read(MPI_Comm comm,
const char *fname);
333 double InnerProduct(HypreParVector &x, HypreParVector &y);
334 double InnerProduct(HypreParVector *x, HypreParVector *y);
339 double ParNormlp(
const Vector &vec,
double p, MPI_Comm comm);
347 hypre_ParCSRMatrix *A;
374 signed char diagOwner, offdOwner, colMapOwner;
377 signed char ParCSROwner;
390 void Write(
MemoryClass mc,
bool set_diag =
true,
bool set_offd =
true);
398 hypre_CSRMatrix *hypre_csr,
404 static signed char CopyBoolCSR(
Table *bool_csr,
406 hypre_CSRMatrix *hypre_csr);
413 static signed char HypreCsrToMem(hypre_CSRMatrix *h_mat,
MemoryType h_mat_mt,
467 bool own_diag_offd =
false);
482 HYPRE_Int *diag_i, HYPRE_Int *diag_j,
double *diag_data,
483 HYPRE_Int *offd_i, HYPRE_Int *offd_j,
double *offd_data,
484 HYPRE_Int offd_num_cols,
486 bool hypre_arrays =
false);
511 HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd,
513 HYPRE_Int cmap_size);
537 operator hypre_ParCSRMatrix*()
const {
return A; }
538 #ifndef HYPRE_PAR_CSR_MATRIX_STRUCT 540 operator HYPRE_ParCSRMatrix() {
return (HYPRE_ParCSRMatrix) A; }
546 void SetOwnerFlags(
signed char diag,
signed char offd,
signed char colmap);
607 bool interleaved_rows =
false,
608 bool interleaved_cols =
false)
const;
615 #if MFEM_HYPRE_VERSION >= 21800 617 double threshold=0.0)
const;
623 return internal::to_int(
624 hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(A)));
630 return internal::to_int(
631 hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(A)));
636 {
return hypre_ParCSRMatrixGlobalNumRows(A); }
640 {
return hypre_ParCSRMatrixGlobalNumCols(A); }
674 double alpha = 1.0,
double beta = 0.0)
const;
676 HYPRE_Int
Mult(HYPRE_ParVector x, HYPRE_ParVector y,
677 double alpha = 1.0,
double beta = 0.0)
const;
684 double alpha = 1.0,
double beta = 0.0)
const;
695 {
Mult(1.0, x, 0.0, y); }
705 {
Mult(
a, x, 1.0, y); }
707 const double a = 1.0)
const override 726 internal::hypre_ParCSRMatrixBooleanMatvec(A,
alpha, const_cast<int*>(x),
736 internal::hypre_ParCSRMatrixBooleanMatvecT(A,
alpha, const_cast<int*>(x),
744 #if MFEM_HYPRE_VERSION < 22200 745 internal::hypre_ParCSRMatrixSetConstantValues(A, value);
747 hypre_ParCSRMatrixSetConstantValues(A, value);
765 MFEM_VERIFY(internal::hypre_ParCSRMatrixSum(A, beta, B.A) == 0,
766 "error in hypre_ParCSRMatrixSum");
879 void Print(
const char *fname, HYPRE_Int offi = 0, HYPRE_Int offj = 0)
const;
881 void Read(MPI_Comm comm,
const char *fname);
913 #if MFEM_HYPRE_VERSION >= 21800 926 const Vector *
b, HypreParVector *d,
933 HypreParMatrix *
Add(
double alpha,
const HypreParMatrix &A,
934 double beta,
const HypreParMatrix &B);
938 HypreParMatrix *
ParMult(
const HypreParMatrix *A,
const HypreParMatrix *B,
939 bool own_matrix =
false);
943 HypreParMatrix *
ParAdd(
const HypreParMatrix *A,
const HypreParMatrix *B);
946 HypreParMatrix *
RAP(
const HypreParMatrix *A,
const HypreParMatrix *P);
948 HypreParMatrix *
RAP(
const HypreParMatrix * Rt,
const HypreParMatrix *A,
949 const HypreParMatrix *P);
960 Array2D<double> *blockCoeff=NULL);
966 void EliminateBC(
const HypreParMatrix &A,
const HypreParMatrix &Ae,
967 const Array<int> &ess_dof_list,
const Vector &X, Vector &B);
1042 #if !defined(HYPRE_USING_GPU) 1140 virtual operator HYPRE_Solver()
const = 0;
1143 virtual HYPRE_PtrToParSolverFcn
SetupFcn()
const = 0;
1145 virtual HYPRE_PtrToParSolverFcn
SolveFcn()
const = 0;
1159 {
mfem_error(
"HypreSolvers do not support SetOperator!"); }
1188 #if MFEM_HYPRE_VERSION >= 21800 1197 virtual operator HYPRE_Solver()
const {
return NULL; }
1200 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSROnProcTriSetup; }
1202 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSROnProcTriSolve; }
1218 HYPRE_Solver pcg_solver;
1252 HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_it);
1253 num_iterations = internal::to_int(num_it);
1258 HYPRE_ParCSRPCGGetFinalRelativeResidualNorm(pcg_solver,
1263 virtual operator HYPRE_Solver()
const {
return pcg_solver; }
1267 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSetup; }
1270 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSolve; }
1283 HYPRE_Solver gmres_solver;
1288 void SetDefaultOptions();
1316 HYPRE_ParCSRGMRESGetNumIterations(gmres_solver, &num_it);
1317 num_iterations = internal::to_int(num_it);
1322 HYPRE_ParCSRGMRESGetFinalRelativeResidualNorm(gmres_solver,
1327 virtual operator HYPRE_Solver()
const {
return gmres_solver; }
1331 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSetup; }
1334 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSolve; }
1347 HYPRE_Solver fgmres_solver;
1352 void SetDefaultOptions();
1379 HYPRE_ParCSRFlexGMRESGetNumIterations(fgmres_solver, &num_it);
1380 num_iterations = internal::to_int(num_it);
1385 HYPRE_ParCSRFlexGMRESGetFinalRelativeResidualNorm(fgmres_solver,
1390 virtual operator HYPRE_Solver()
const {
return fgmres_solver; }
1394 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRFlexGMRESSetup; }
1397 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRFlexGMRESSolve; }
1410 virtual operator HYPRE_Solver()
const {
return NULL; }
1413 {
return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentitySetup; }
1415 {
return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentity; }
1426 virtual operator HYPRE_Solver()
const {
return NULL; }
1431 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScaleSetup; }
1433 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScale; }
1448 HYPRE_Solver sai_precond;
1451 void SetDefaultOptions();
1456 void ResetSAIPrecond(MPI_Comm comm);
1465 void SetParams(
double threshold,
int max_levels);
1473 virtual operator HYPRE_Solver()
const {
return sai_precond; }
1476 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSetup; }
1478 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSolve; }
1494 HYPRE_Solver euc_precond;
1497 void SetDefaultOptions();
1502 void ResetEuclidPrecond(MPI_Comm comm);
1518 virtual operator HYPRE_Solver()
const {
return euc_precond; }
1521 {
return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSetup; }
1523 {
return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSolve; }
1528 #if MFEM_HYPRE_VERSION >= 21900 1546 HYPRE_Solver ilu_precond;
1549 void SetDefaultOptions();
1555 void ResetILUPrecond();
1566 void SetType(HYPRE_Int ilu_type);
1568 void SetTol(HYPRE_Real tol);
1575 virtual operator HYPRE_Solver()
const {
return ilu_precond; }
1581 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ILUSetup; }
1585 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ILUSolve; }
1593 HYPRE_Solver amg_precond;
1602 void RecomputeRBMs();
1605 void SetDefaultOptions();
1610 void ResetAMGPrecond();
1630 #if MFEM_HYPRE_VERSION >= 21800 1640 const std::string &postrelax=
"FFC");
1644 { HYPRE_BoomerAMGSetStrongThresholdR(amg_precond, strengthR); }
1648 { HYPRE_BoomerAMGSetFilterThresholdR(amg_precond, filterR); }
1652 { HYPRE_BoomerAMGSetRestriction(amg_precond, restrict_type); }
1656 { HYPRE_BoomerAMGSetIsTriangular(amg_precond, 1); }
1660 { HYPRE_BoomerAMGSetGMRESSwitchR(amg_precond, gmres_switch); }
1665 HYPRE_BoomerAMGSetCycleNumSweeps(amg_precond, prerelax, 1);
1666 HYPRE_BoomerAMGSetCycleNumSweeps(amg_precond, postrelax, 2);
1671 { HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level); }
1674 { HYPRE_BoomerAMGSetMaxIter(amg_precond, max_iter); }
1678 { HYPRE_BoomerAMGSetMaxLevels(amg_precond, max_levels); }
1682 { HYPRE_BoomerAMGSetTol(amg_precond, tol); }
1686 { HYPRE_BoomerAMGSetStrongThreshold(amg_precond, strength); }
1690 { HYPRE_BoomerAMGSetInterpType(amg_precond, interp_type); }
1694 { HYPRE_BoomerAMGSetCoarsenType(amg_precond, coarsen_type); }
1698 { HYPRE_BoomerAMGSetRelaxType(amg_precond, relax_type); }
1702 { HYPRE_BoomerAMGSetCycleType(amg_precond, cycle_type); }
1707 HYPRE_BoomerAMGGetNumIterations(amg_precond, &num_it);
1708 num_iterations = internal::to_int(num_it);
1714 HYPRE_BoomerAMGSetNumFunctions(amg_precond, blocksize);
1715 HYPRE_BoomerAMGSetNodal(amg_precond, 1);
1720 { HYPRE_BoomerAMGSetAggNumLevels(amg_precond, num_levels); }
1723 virtual operator HYPRE_Solver()
const {
return amg_precond; }
1726 {
return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSetup; }
1728 {
return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSolve; }
1736 HypreParMatrix*
DiscreteGrad(ParFiniteElementSpace *edge_fespace,
1737 ParFiniteElementSpace *vert_fespace);
1739 HypreParMatrix*
DiscreteCurl(ParFiniteElementSpace *face_fespace,
1740 ParFiniteElementSpace *edge_fespace);
1751 void MakeSolver(
int sdim,
int cycle_type);
1768 int ams_cycle_type = 0;
1772 bool singular =
false;
1774 int print_level = 1;
1778 void ResetAMSPrecond();
1804 HYPRE_AMSSetBetaPoissonMatrix(ams, NULL);
1809 virtual operator HYPRE_Solver()
const {
return ams; }
1812 {
return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSetup; }
1814 {
return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSolve; }
1848 const int cycle_type = 11;
1850 const int ams_cycle_type = 14;
1852 int print_level = 1;
1856 void ResetADSPrecond();
1876 virtual operator HYPRE_Solver()
const {
return ads; }
1879 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSetup; }
1881 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSolve; }
1921 HYPRE_Solver lobpcg_solver;
1924 mv_InterfaceInterpreter interpreter;
1927 HYPRE_MatvecFunctions matvec_fn;
1933 class HypreMultiVector;
1936 HypreMultiVector * multi_vec;
1945 class HypreMultiVector
1949 mv_MultiVectorPtr mv_ptr;
1959 mv_InterfaceInterpreter & interpreter);
1960 ~HypreMultiVector();
1963 void Randomize(HYPRE_Int seed);
1971 operator mv_MultiVectorPtr()
const {
return mv_ptr; }
1973 mv_MultiVectorPtr & GetMultiVector() {
return mv_ptr; }
1976 static void * OperatorMatvecCreate(
void *A,
void *x );
1977 static HYPRE_Int OperatorMatvec(
void *matvec_data,
1978 HYPRE_Complex
alpha,
1983 static HYPRE_Int OperatorMatvecDestroy(
void *matvec_data );
1985 static HYPRE_Int PrecondSolve(
void *solver,
1989 static HYPRE_Int PrecondSetup(
void *solver,
2057 HYPRE_Solver ame_solver;
2063 HYPRE_Real * eigenvalues;
2066 HYPRE_ParVector * multi_vec;
2071 void createDummyVectors()
const;
2103 #endif // MFEM_USE_MPI
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
MemoryType GetHypreMemoryType()
The MemoryType used by MFEM when allocating arrays for Hypre objects.
virtual ~HypreBoomerAMG()
void SetAbsTol(double atol)
void GetEigenvalues(Array< double > &eigenvalues) const
Collect the converged eigenvalues.
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
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 ...
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
void SetType(HYPRE_Int ilu_type)
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's GMRES.
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
HypreADS(ParFiniteElementSpace *face_fespace)
HypreParMatrix * ExtractSubmatrix(const Array< int > &indices, double threshold=0.0) const
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(), hypre will be finalized automatically at program exit.
static constexpr Type default_type
HypreEuclid(MPI_Comm comm)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void EliminateBC(const HypreParMatrix &A, const HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
Eliminate essential BC specified by ess_dof_list from the solution X to the r.h.s. B.
The Auxiliary-space Maxwell Solver in hypre.
HypreParVector * X0
FIR Filter Temporary Vectors.
void SetStrengthThresh(double strength)
Expert option - consult hypre documentation/team.
double min_eig_est
Minimal eigenvalue estimate for polynomial smoothing.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
static MemoryClass GetHostMemoryClass()
Get the current Host MemoryClass. This is the MemoryClass used by most MFEM host Memory objects...
const Memory< double > & GetDiagMemoryData() const
ErrorMode
How to treat errors returned by hypre function calls.
The Auxiliary-space Divergence Solver in hypre.
void SetPrintLevel(int print_lvl)
MFEM_DEPRECATED HypreParMatrix * GetData()
Deprecated. Use HypreTriSolve::GetData() const instead.
void SetRowScale(int row_scale)
void SetMassMatrix(const HypreParMatrix &M)
void WrapMemoryWrite(Memory< double > &mem)
Replace the HypreParVector's data with the given Memory, mem, and prepare the vector for write access...
int setup_called
Was hypre's Setup function called already?
Device memory; using CUDA or HIP *Malloc and *Free.
const HypreParMatrix * GetData() const
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
void SetLogging(int logging)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
const HypreParVector & GetEigenvector(unsigned int i) const
Extract a single eigenvector.
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to 'master'.
void Read_IJMatrix(MPI_Comm comm, const char *fname)
Read a matrix saved as a HYPRE_IJMatrix.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
GMRES Solve function.
static void Finalize()
Finalize hypre (called automatically at program exit if Hypre::Init() has been called).
HypreParMatrix & Add(const double beta, const HypreParMatrix &B)
HypreParVector * B
Right-hand side and solution vector.
Wrapper for Hypre's native parallel ILU preconditioner.
double window_params[3]
Parameters for windowing function of FIR filter.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
PCG Setup function.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
FGMRES Solve function.
void GetEigenvalues(Array< double > &eigenvalues) const
Collect the converged eigenvalues.
void SetInitialVectors(int num_vecs, HypreParVector **vecs)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void SetWindowByName(const char *window_name)
Convenience function for setting canonical windowing parameters.
HypreTriSolve(const HypreParMatrix &A)
Issue warnings on hypre errors.
T * Write(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for write access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true...
void SetPreconditioner(HypreSolver &precond)
HypreILU()
Constructor; sets the default options.
void GetNumIterations(int &num_iterations) const
int GetNumCols() const
Returns the number of columns in the diagonal block of the ParCSRMatrix.
void SetRestriction(int restrict_type)
Expert option - consult hypre documentation/team.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
void AssembleDiagonal(Vector &diag) const override
Return the diagonal of the matrix (Operator interface).
HypreAMS(ParFiniteElementSpace *edge_fespace)
Construct the AMS solver on the given edge finite element space.
void SetCycleNumSweeps(int prerelax, int postrelax)
Expert option - consult hypre documentation/team.
void SetPrintLevel(int print_lvl)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void SetPreconditioner(Solver &precond)
MFEM_DEPRECATED void SetZeroInintialIterate()
deprecated: use SetZeroInitialIterate()
void ResetTranspose() const
Reset (destroy) the internal transpose matrix that is created by EnsureMultTranspose() and MultTransp...
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Memory< HYPRE_Int > & GetDiagMemoryJ()
void SetMassMatrix(Operator &M)
Abstract parallel finite element space.
HypreParVector()
Default constructor, no underlying hypre_ParVector is created.
const Memory< HYPRE_Int > & GetDiagMemoryI() const
void SetPrintLevel(int logging)
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void HostRead() const
Update the internal hypre_ParCSRMatrix object, A, to be on host.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
void SetAdvectiveOptions(int distance=15, const std::string &prerelax="", const std::string &postrelax="FFC")
HYPRE_BigInt N() const
Returns the global number of columns.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
FGMRES Setup function.
void BlockInverseScale(const HypreParMatrix *A, HypreParMatrix *C, const Vector *b, HypreParVector *d, int blocksize, BlockInverseScaleJob job)
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
void HostWrite()
Update the internal hypre_ParCSRMatrix object, A, to be on host.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
HypreParMatrix & operator+=(const HypreParMatrix &B)
bool pos_l1_norms
If set, take absolute values of the computed l1_norms.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
HYPRE_BigInt NNZ() const
Returns the global number of nonzeros.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
void SetPrintLevel(int print_lvl)
void SetZeroInitialIterate()
non-hypre setting
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
void SetStrongThresholdR(double strengthR)
Expert option - consult hypre documentation/team.
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
void SetParams(double threshold, int max_levels)
double poly_fraction
Fraction of spectrum to smooth for polynomial relaxation.
int GetOwnership() const
Gets ownership of the internal hypre_ParVector.
void SetPrintLevel(int print_lvl)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
HypreParMatrix * DiscreteGrad(ParFiniteElementSpace *edge_fespace, ParFiniteElementSpace *vert_fespace)
Compute the discrete gradient matrix between the nodal linear and ND1 spaces.
int poly_scale
Apply the polynomial smoother to A or D^{-1/2} A D^{-1/2}.
void MergeDiagAndOffd(SparseMatrix &merged)
Get a single SparseMatrix containing all rows from this processor, merged from the diagonal and off-d...
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
GMRES Setup function.
void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
The "Boolean" analog of y = alpha * A^T * x + beta * y, where elements in the sparsity pattern of the...
void SetMaxIter(HYPRE_Int max_iter)
HypreFGMRES(MPI_Comm comm)
The identity operator as a hypre solver.
void SetWindowParameters(double a, double b, double c)
Set parameters for windowing function for FIR smoother.
void EliminateZeroRows()
If a row contains only zeros, set its diagonal to 1.
void SetLoadBal(double loadbal)
HypreGMRES(MPI_Comm comm)
const HYPRE_BigInt * ColPart() const
Returns the column partitioning (const version)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
void HypreWrite()
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
HypreLOBPCG(MPI_Comm comm)
HYPRE_BigInt GlobalSize() const
Returns the global number of rows.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
PCG Solve function.
const HYPRE_BigInt * Partitioning() const
Returns the parallel row/column partitioning.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
MFEM_DEPRECATED HYPRE_BigInt * Partitioning()
Returns a non-const pointer to the parallel row/column partitioning. Deprecated in favor of HypreParV...
void EliminateBC(const HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B) const
Eliminate essential BC specified by ess_dof_list from the solution X to the r.h.s. B.
hypre_ParVector * StealParVector()
Changes the ownership of the vector.
The BoomerAMG solver in hypre.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void SetRelaxType(int relax_type)
Expert option - consult hypre documentation/team.
void SetSymmetry(int sym)
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
void HypreStealOwnership(HypreParMatrix &A_hyp, SparseMatrix &A_diag)
Make A_hyp steal ownership of its diagonal part A_diag.
virtual ~HypreParaSails()
The ParaSails preconditioner in hypre.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
void SetMaxIter(int max_iter)
double * l1_norms
l1 norms of the rows of A
void SetOwnerFlags(signed char diag, signed char offd, signed char colmap)
Explicitly set the three ownership flags, see docs for diagOwner etc.
virtual void Setup(const HypreParVector &b, HypreParVector &x) const
Set up the solver (if not set up already, also called automatically by HypreSolver::Mult).
void SetLogging(int logging)
Jacobi preconditioner in hypre.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const =0
hypre's internal Solve function
void SetRelTol(double rel_tol)
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0)
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
int to_int(const std::string &str)
Convert a string to an int.
Abort on hypre errors (default in base class)
void SetLogging(int logging)
void SetPrintLevel(int print_level)
void Print(const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0) const
Prints the locally owned rows in parallel.
void SetOperatorSymmetry(bool is_sym)
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
void GetFinalResidualNorm(double &final_res_norm) const
HypreParaSails(MPI_Comm comm)
void HypreWrite()
Prepare the HypreParVector for write access in hypre's device memory space, HYPRE_MEMORY_DEVICE.
void WrapHypreParVector(hypre_ParVector *y, bool owner=true)
Converts hypre's format to HypreParVector.
void AbsMultTranspose(double a, const Vector &x, double b, Vector &y) const
Computes y = a * |At| * x + b * y, using entry-wise absolute values of the transpose of the matrix A...
void SetTol(double tol)
Expert option - consult hypre documentation/team.
bool WrapVectors(const Vector &b, Vector &x) const
Makes the internal HypreParVectors B and X wrap the input vectors b and x.
void Read(MPI_Comm comm, const char *fname)
Reads a HypreParVector from files saved with HypreParVector::Print.
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_BigInt *row_starts=NULL) const
Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result a...
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
const HypreParMatrix * GetData() const
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
const HypreParMatrix * A
The linear system matrix.
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.
void SetFilterThresholdR(double filterR)
Expert option - consult hypre documentation/team.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's FGMRES.
void GetFinalResidualNorm(double &final_res_norm) const
virtual void MultTranspose(const Vector &b, Vector &x) const
Apply transpose of the smoother to relax the linear system Ax=b.
Ignore hypre errors (see e.g. HypreADS)
double relax_weight
Damping coefficient (usually <= 1)
HYPRE_BigInt * GetColStarts() const
Return the parallel column partitioning array.
Parallel smoothers in hypre.
void SetElasticityOptions(ParFiniteElementSpace *fespace)
void AddMult(const Vector &x, Vector &y, const double a=1.0) const override
Operator application: y+=A(x) (default) or y+=a*A(x).
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
const T * Read(const Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true...
void SetZeroInitialIterate()
non-hypre setting
void SetLevelOfFill(HYPRE_Int lev_fill)
Set the fill level for ILU(k); the default is k=1.
HypreParMatrix * A
The linear system matrix.
const Memory< HYPRE_Int > & GetDiagMemoryJ() const
void SetCoarsening(int coarsen_type)
Expert option - consult hypre documentation/team.
void SetGMRESSwitchR(int gmres_switch)
Expert option - consult hypre documentation/team.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
ILU Solve function.
void AddMultTranspose(const Vector &x, Vector &y, const double a=1.0) const override
Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x).
Dynamic 2D array using row-major layout.
virtual void SetOperator(const Operator &op)
virtual HYPRE_PtrToParSolverFcn SetupFcn() const =0
hypre's internal Setup function
void SetLogging(int logging)
HypreParMatrix()
An empty matrix to be used as a reference to an existing matrix.
bool A_is_symmetric
A flag that indicates whether the linear system matrix A is symmetric.
void EnsureMultTranspose() const
Ensure the action of the transpose is performed fast.
void SetMaxIter(int max_iter)
void SetMaxIter(int max_iter)
void SetMaxIter(int max_iter)
void SetPolyOptions(int poly_order, double poly_fraction, int eig_est_cg_iter=10)
Set parameters for polynomial smoothing.
HypreDiagScale(const HypreParMatrix &A)
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
MPI_Comm GetComm() const
MPI communicator.
void BooleanMult(int alpha, const int *x, int beta, int *y)
The "Boolean" analog of y = alpha * A * x + beta * y, where elements in the sparsity pattern of the m...
void SetInterpolation(int interp_type)
Expert option - consult hypre documentation/team.
Wrapper for hypre's parallel vector class.
HypreParMatrix * EliminateCols(const Array< int > &cols)
Type
Enumeration defining IDs for some classes derived from Operator.
double p(const Vector &x, double t)
static MemoryType GetHostMemoryType()
Get the current Host MemoryType. This is the MemoryType used by most MFEM classes when allocating mem...
void PrintHash(std::ostream &out) const
Print sizes and hashes for all data arrays of the HypreParMatrix from the local MPI rank...
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
MemoryType
Memory types supported by MFEM.
void HypreRead() const
Prepare the HypreParVector for read access in hypre's device memory space, HYPRE_MEMORY_DEVICE.
void SetPrintLevel(int print_lvl)
void PrintCommPkg(std::ostream &out=mfem::out) const
Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
signed char OwnsColMap() const
Get colmap ownership flag.
void Solve()
Solve the eigenproblem.
void SetData(double *data_)
Sets the data of the Vector and the hypre_ParVector to data_.
const HYPRE_BigInt * RowPart() const
Returns the row partitioning (const version)
void SetNumModes(int num_eigs)
HypreParMatrix(hypre_ParCSRMatrix *a, bool owner=true)
Converts hypre's format to HypreParMatrix.
virtual ~HypreDiagScale()
void MultTranspose(const Vector &x, Vector &y) const override
Computes y = A^t * x.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A + B.
double max_eig_est
Maximal eigenvalue estimate for polynomial smoothing.
A class to initialize the size of a Tensor.
void SetRelTol(double rel_tol)
MPI_Comm GetComm() const
MPI communicator.
HypreParMatrix * Transpose() const
Returns the transpose of *this.
void GetBlocks(Array2D< HypreParMatrix *> &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
void SetMaxLevels(int max_levels)
Expert option - consult hypre documentation/team.
void SetMaxIter(int max_iter)
void SetFIRCoefficients(double max_eig)
Compute window and Chebyshev coefficients for given polynomial order.
void SetRandomSeed(int s)
void EliminateRows(const Array< int > &rows)
Eliminate rows from the diagonal and off-diagonal blocks of the matrix.
double * fir_coeffs
Combined coefficients for windowing and Chebyshev polynomials.
void SetIsTriangular()
Expert option - consult hypre documentation/team.
void Threshold(double threshold=0.0)
Remove values smaller in absolute value than some threshold.
double InnerProduct(HypreParVector *x, HypreParVector *y)
void SetErrorMode(ErrorMode err_mode) const
Set the behavior for treating hypre errors, see the ErrorMode enum. The default mode in the base clas...
T * ReadWrite(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read+write access to mem with the mfem::Device's DeviceMemoryClass, if on_dev = true, or the mfem::Device's HostMemoryClass, otherwise.
void SetPrintLevel(HYPRE_Int print_level)
Set the print level: 0 = none, 1 = setup, 2 = solve, 3 = setup+solve.
void GetNumIterations(int &num_iterations) const
Host memory; using new[] and delete[].
void SetSubSpaceProjector(Operator &proj)
constexpr MemoryClass GetHypreMemoryClass()
The MemoryClass used by Hypre objects.
void Print(const char *fname) const
Prints the locally owned rows in parallel.
HYPRE_BigInt * ColPart()
Returns the column partitioning.
void HypreReadWrite()
Prepare the HypreParVector for read and write access in hypre's device memory space, HYPRE_MEMORY_DEVICE.
Memory< HYPRE_Int > & GetDiagMemoryI()
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
MFEM_DEPRECATED void SetZeroInintialIterate()
deprecated: use SetZeroInitialIterate()
void SetTaubinOptions(double lambda, double mu, int iter)
Set parameters for Taubin's lambda-mu method.
HypreParVector * B
Right-hand side and solution vectors.
HYPRE_BigInt * RowPart()
Returns the row partitioning.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
A simple singleton class for hypre's global settings, that 1) calls HYPRE_Init() and sets some GPU-re...
DiagonalPolicy
Defines operator diagonal policy upon elimination of rows and/or columns.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
void GetNumIterations(int &num_iterations) const
void SetZeroInitialIterate()
non-hypre setting
void SetOperator(Operator &A)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
void SetLocalReordering(HYPRE_Int reorder_type)
void HypreRead() const
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
int relax_times
Number of relaxation sweeps.
void SetCycleType(int cycle_type)
Expert option - consult hypre documentation/team.
void operator*=(double s)
Scale all entries by s: A_scaled = s*A.
void SetPrintLevel(int logging)
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...
int eig_est_cg_iter
Number of CG iterations to determine eigenvalue estimates.
void SetOwnership(int own)
Sets ownership of the internal hypre_ParVector.
Abstract class for hypre's solvers and preconditioners.
Flexible GMRES solver in hypre.
ErrorMode error_mode
How to treat hypre errors.
void AbsMult(double a, const Vector &x, double b, Vector &y) const
Computes y = a * |A| * x + b * y, using entry-wise absolute values of the matrix A.
HypreParVector & operator=(double d)
Set constant values.
void GetNumIterations(int &num_iterations) const
void SetAbsTol(double tol)
signed char OwnsDiag() const
Get diag ownership flag.
HypreParMatrix & operator=(double value)
Initialize all entries with value.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
ILU Setup function.
hypre_ParCSRMatrix * StealData()
Changes the ownership of the matrix.
void SetPrecondUsageMode(int pcg_mode)
void SetOperator(const HypreParMatrix &A)
void SetSystemsOptions(int dim, bool order_bynodes=false)
virtual ~HypreParMatrix()
Calls hypre's destroy function.
ID for class HypreParMatrix.
void SetFilter(double filter)
void DropSmallEntries(double tol)
Wrapper for hypre_ParCSRMatrixDropSmallEntries in different versions of hypre. Drop off-diagonal entr...
void SetAggressiveCoarsening(int num_levels)
Expert option - consult hypre documentation/team.
void SetMaxIter(int max_iter)
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
Vector * GlobalVector() const
Returns the global vector in each processor.
void GetOffd(SparseMatrix &offd, HYPRE_BigInt *&cmap) const
Get the local off-diagonal block. NOTE: 'offd' will not own any data.
double omega
SOR parameter (usually in (0,2))
HYPRE_BigInt M() const
Returns the global number of rows.
signed char OwnsOffd() const
Get offd ownership flag.
MemoryClass GetMemoryClass() const override
Return the MemoryClass preferred by the Operator.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A^t * x + beta * y.
void SetNodal(int blocksize)
Expert option - consult hypre documentation/team.
void WrapHypreParCSRMatrix(hypre_ParCSRMatrix *a, bool owner=true)
Converts hypre's format to HypreParMatrix.
void GetFinalResidualNorm(double &final_res_norm) const
Wrapper for hypre's ParCSR matrix class.
void HostReadWrite()
Update the internal hypre_ParCSRMatrix object, A, to be on host.
Memory< double > auxB
Auxiliary buffers for the case when the input or output arrays in methods like Mult(const Vector &...
~HypreParVector()
Calls hypre's destroy function.
void SetNumModes(int num_eigs)
HypreParMatrix * HypreParMatrixFromBlocks(Array2D< HypreParMatrix *> &blocks, Array2D< double > *blockCoeff)
Returns a merged hypre matrix constructed from hypre matrix blocks.
Memory< double > & GetDiagMemoryData()
void SetSingularProblem()
Set this option when solving a curl-curl problem with zero mass term.
void Solve()
Solve the eigenproblem.
HypreParVector * V
Temporary vectors.
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
MemoryClass
Memory classes identify sets of memory types.
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
void WrapMemoryReadWrite(Memory< double > &mem)
Replace the HypreParVector's data with the given Memory, mem, and prepare the vector for read and wri...
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A * x + beta * y.
HypreParVector CreateCompatibleVector() const
Constructs a HypreParVector compatible with the calling vector.
int poly_order
Order of the smoothing polynomial.
void SetTol(HYPRE_Real tol)
int GetNumRows() const
Returns the number of rows in the diagonal block of the ParCSRMatrix.
double lambda
Taubin's lambda-mu method parameters.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
void WrapMemoryRead(const Memory< double > &mem)
Replace the HypreParVector's data with the given Memory, mem, and prepare the vector for read access ...
MFEM_DEPRECATED void SetZeroInintialIterate()
deprecated: use SetZeroInitialIterate()
const HypreParVector & GetEigenvector(unsigned int i) const
Extract a single eigenvector.
MFEM_DEPRECATED HypreParMatrix * GetData()
Deprecated. Use HypreDiagScale::GetData() const instead.
void HypreReadWrite()
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
HYPRE_BigInt * GetRowStarts() const
Return the parallel row partitioning array.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function