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_CUDA) && !defined(MFEM_USE_CUDA)
36 #error "MFEM_USE_CUDA=YES is required when HYPRE is built with CUDA!"
45 class ParFiniteElementSpace;
51 template <
typename int_type>
52 inline int to_int(int_type i)
54 MFEM_ASSERT(int_type(
int(i)) == i,
"overflow converting int_type to int");
59 template <>
inline int to_int(
int i) {
return i; }
64 inline int to_int(HYPRE_Int i)
66 MFEM_ASSERT(HYPRE_Int(
int(i)) == i,
"overflow converting HYPRE_Int to int");
77 #ifndef HYPRE_USING_CUDA
79 #elif defined(HYPRE_USING_UNIFIED_MEMORY)
89 #ifndef HYPRE_USING_CUDA
91 #elif defined(HYPRE_USING_UNIFIED_MEMORY)
110 inline void _SetDataAndSize_();
117 own_ParVector =
false;
174 operator hypre_ParVector*()
const {
return x; }
175 #ifndef HYPRE_PAR_VECTOR_STRUCT
177 operator HYPRE_ParVector()
const {
return (HYPRE_ParVector) x; }
253 void Print(
const char *fname)
const;
258 #ifdef MFEM_USE_SUNDIALS
261 MFEM_DEPRECATED
virtual N_Vector
ToNVector();
267 double InnerProduct(HypreParVector &x, HypreParVector &y);
268 double InnerProduct(HypreParVector *x, HypreParVector *y);
273 double ParNormlp(
const Vector &vec,
double p, MPI_Comm comm);
281 hypre_ParCSRMatrix *A;
308 signed char diagOwner, offdOwner, colMapOwner;
311 signed char ParCSROwner;
324 void Write(
MemoryClass mc,
bool set_diag =
true,
bool set_offd =
true);
332 hypre_CSRMatrix *hypre_csr,
338 static signed char CopyBoolCSR(
Table *bool_csr,
340 hypre_CSRMatrix *hypre_csr);
344 static void CopyCSR_J(hypre_CSRMatrix *hypre_csr,
int *J);
351 static signed char HypreCsrToMem(hypre_CSRMatrix *h_mat,
MemoryType h_mat_mt,
405 bool own_diag_offd =
false);
420 HYPRE_Int *diag_i, HYPRE_Int *diag_j,
double *diag_data,
421 HYPRE_Int *offd_i, HYPRE_Int *offd_j,
double *offd_data,
422 HYPRE_Int offd_num_cols,
424 bool hypre_arrays =
false);
449 HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd,
451 HYPRE_Int cmap_size);
475 operator hypre_ParCSRMatrix*()
const {
return A; }
476 #ifndef HYPRE_PAR_CSR_MATRIX_STRUCT
478 operator HYPRE_ParCSRMatrix() {
return (HYPRE_ParCSRMatrix) A; }
484 void SetOwnerFlags(
signed char diag,
signed char offd,
signed char colmap);
545 bool interleaved_rows =
false,
546 bool interleaved_cols =
false)
const;
553 #if MFEM_HYPRE_VERSION >= 21800
555 double threshold=0.0)
const;
562 hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(A)));
569 hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(A)));
574 {
return hypre_ParCSRMatrixGlobalNumRows(A); }
578 {
return hypre_ParCSRMatrixGlobalNumCols(A); }
594 double alpha = 1.0,
double beta = 0.0)
const;
596 HYPRE_Int
Mult(HYPRE_ParVector x, HYPRE_ParVector y,
597 double alpha = 1.0,
double beta = 0.0)
const;
600 double alpha = 1.0,
double beta = 0.0)
const;
606 {
Mult(1.0, x, 0.0, y); }
623 internal::hypre_ParCSRMatrixBooleanMatvec(A, alpha, const_cast<int*>(x),
633 internal::hypre_ParCSRMatrixBooleanMatvecT(A, alpha, const_cast<int*>(x),
641 #if MFEM_HYPRE_VERSION < 22200
642 internal::hypre_ParCSRMatrixSetConstantValues(A, value);
644 hypre_ParCSRMatrixSetConstantValues(A, value);
662 MFEM_VERIFY(internal::hypre_ParCSRMatrixSum(A, beta, B.A) == 0,
663 "error in hypre_ParCSRMatrixSum");
760 void Print(
const char *fname, HYPRE_Int offi = 0, HYPRE_Int offj = 0)
const;
762 void Read(MPI_Comm comm,
const char *fname);
782 #if MFEM_HYPRE_VERSION >= 21800
795 const Vector *
b, HypreParVector *d,
802 HypreParMatrix *
Add(
double alpha,
const HypreParMatrix &A,
803 double beta,
const HypreParMatrix &B);
807 HypreParMatrix *
ParMult(
const HypreParMatrix *A,
const HypreParMatrix *B,
808 bool own_matrix =
false);
812 HypreParMatrix *
ParAdd(
const HypreParMatrix *A,
const HypreParMatrix *B);
815 HypreParMatrix *
RAP(
const HypreParMatrix *A,
const HypreParMatrix *P);
817 HypreParMatrix *
RAP(
const HypreParMatrix * Rt,
const HypreParMatrix *A,
818 const HypreParMatrix *P);
829 Array2D<double> *blockCoeff=NULL);
835 void EliminateBC(
const HypreParMatrix &A,
const HypreParMatrix &Ae,
836 const Array<int> &ess_dof_list,
const Vector &X, Vector &B);
911 #ifndef HYPRE_USING_CUDA
1002 virtual operator HYPRE_Solver()
const = 0;
1005 virtual HYPRE_PtrToParSolverFcn
SetupFcn()
const = 0;
1007 virtual HYPRE_PtrToParSolverFcn
SolveFcn()
const = 0;
1010 {
mfem_error(
"HypreSolvers do not support SetOperator!"); }
1033 #if MFEM_HYPRE_VERSION >= 21800
1042 virtual operator HYPRE_Solver()
const {
return NULL; }
1045 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSROnProcTriSetup; }
1047 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSROnProcTriSolve; }
1063 HYPRE_Solver pcg_solver;
1097 HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_it);
1102 virtual operator HYPRE_Solver()
const {
return pcg_solver; }
1106 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSetup; }
1109 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSolve; }
1122 HYPRE_Solver gmres_solver;
1127 void SetDefaultOptions();
1153 virtual operator HYPRE_Solver()
const {
return gmres_solver; }
1157 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSetup; }
1160 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSolve; }
1173 HYPRE_Solver fgmres_solver;
1178 void SetDefaultOptions();
1203 virtual operator HYPRE_Solver()
const {
return fgmres_solver; }
1207 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRFlexGMRESSetup; }
1210 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRFlexGMRESSolve; }
1223 virtual operator HYPRE_Solver()
const {
return NULL; }
1226 {
return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentitySetup; }
1228 {
return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentity; }
1239 virtual operator HYPRE_Solver()
const {
return NULL; }
1244 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScaleSetup; }
1246 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScale; }
1261 HYPRE_Solver sai_precond;
1264 void SetDefaultOptions();
1269 void ResetSAIPrecond(MPI_Comm comm);
1281 virtual operator HYPRE_Solver()
const {
return sai_precond; }
1284 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSetup; }
1286 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSolve; }
1302 HYPRE_Solver euc_precond;
1305 void SetDefaultOptions();
1310 void ResetEuclidPrecond(MPI_Comm comm);
1320 virtual operator HYPRE_Solver()
const {
return euc_precond; }
1323 {
return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSetup; }
1325 {
return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSolve; }
1330 #if MFEM_HYPRE_VERSION >= 21900
1348 HYPRE_Solver ilu_precond;
1351 void SetDefaultOptions();
1357 void ResetILUPrecond();
1372 virtual operator HYPRE_Solver()
const {
return ilu_precond; }
1378 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ILUSetup; }
1382 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ILUSolve; }
1390 HYPRE_Solver amg_precond;
1399 void RecomputeRBMs();
1402 void SetDefaultOptions();
1407 void ResetAMGPrecond();
1427 #if MFEM_HYPRE_VERSION >= 21800
1437 const std::string &postrelax=
"FFC");
1441 { HYPRE_BoomerAMGSetStrongThresholdR(amg_precond, strengthR); }
1445 { HYPRE_BoomerAMGSetFilterThresholdR(amg_precond, filterR); }
1449 { HYPRE_BoomerAMGSetRestriction(amg_precond, restrict_type); }
1453 { HYPRE_BoomerAMGSetIsTriangular(amg_precond, 1); }
1457 { HYPRE_BoomerAMGSetGMRESSwitchR(amg_precond, gmres_switch); }
1462 HYPRE_BoomerAMGSetCycleNumSweeps(amg_precond, prerelax, 1);
1463 HYPRE_BoomerAMGSetCycleNumSweeps(amg_precond, postrelax, 2);
1468 { HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level); }
1471 { HYPRE_BoomerAMGSetMaxIter(amg_precond, max_iter); }
1475 { HYPRE_BoomerAMGSetMaxLevels(amg_precond, max_levels); }
1479 { HYPRE_BoomerAMGSetTol(amg_precond, tol); }
1483 { HYPRE_BoomerAMGSetStrongThreshold(amg_precond, strength); }
1487 { HYPRE_BoomerAMGSetInterpType(amg_precond, interp_type); }
1491 { HYPRE_BoomerAMGSetCoarsenType(amg_precond, coarsen_type); }
1495 { HYPRE_BoomerAMGSetRelaxType(amg_precond, relax_type); }
1499 { HYPRE_BoomerAMGSetCycleType(amg_precond, cycle_type); }
1504 HYPRE_BoomerAMGGetNumIterations(amg_precond, &num_it);
1511 HYPRE_BoomerAMGSetNumFunctions(amg_precond, blocksize);
1512 HYPRE_BoomerAMGSetNodal(amg_precond, 1);
1517 { HYPRE_BoomerAMGSetAggNumLevels(amg_precond, num_levels); }
1520 virtual operator HYPRE_Solver()
const {
return amg_precond; }
1523 {
return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSetup; }
1525 {
return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSolve; }
1533 HypreParMatrix*
DiscreteGrad(ParFiniteElementSpace *edge_fespace,
1534 ParFiniteElementSpace *vert_fespace);
1536 HypreParMatrix*
DiscreteCurl(ParFiniteElementSpace *face_fespace,
1537 ParFiniteElementSpace *edge_fespace);
1568 virtual operator HYPRE_Solver()
const {
return ams; }
1571 {
return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSetup; }
1573 {
return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSolve; }
1608 virtual operator HYPRE_Solver()
const {
return ads; }
1611 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSetup; }
1613 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSolve; }
1653 HYPRE_Solver lobpcg_solver;
1656 mv_InterfaceInterpreter interpreter;
1659 HYPRE_MatvecFunctions matvec_fn;
1665 class HypreMultiVector;
1668 HypreMultiVector * multi_vec;
1677 class HypreMultiVector
1681 mv_MultiVectorPtr mv_ptr;
1691 mv_InterfaceInterpreter & interpreter);
1692 ~HypreMultiVector();
1695 void Randomize(HYPRE_Int seed);
1703 operator mv_MultiVectorPtr()
const {
return mv_ptr; }
1705 mv_MultiVectorPtr & GetMultiVector() {
return mv_ptr; }
1708 static void * OperatorMatvecCreate(
void *A,
void *x );
1709 static HYPRE_Int OperatorMatvec(
void *matvec_data,
1710 HYPRE_Complex
alpha,
1715 static HYPRE_Int OperatorMatvecDestroy(
void *matvec_data );
1717 static HYPRE_Int PrecondSolve(
void *solver,
1721 static HYPRE_Int PrecondSetup(
void *solver,
1789 HYPRE_Solver ame_solver;
1795 HYPRE_Real * eigenvalues;
1798 HYPRE_ParVector * multi_vec;
1803 void createDummyVectors()
const;
1835 #endif // MFEM_USE_MPI
MemoryType GetHypreMemoryType()
The MemoryType used by MFEM when allocating arrays for Hypre objects.
virtual ~HypreBoomerAMG()
void SetAbsTol(double atol)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
HypreADS(ParFiniteElementSpace *face_fespace)
const HypreParMatrix * GetData() const
static constexpr Type default_type
MPI_Comm GetComm() const
MPI communicator.
HypreEuclid(MPI_Comm comm)
void SetErrorMode(ErrorMode err_mode) const
Set the behavior for treating hypre errors, see the ErrorMode enum. The default mode in the base clas...
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Vector * GlobalVector() const
Returns the global vector in each processor.
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...
ErrorMode
How to treat errors returned by hypre function calls.
The Auxiliary-space Divergence Solver in hypre.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
FGMRES Solve function.
void SetPrintLevel(int print_lvl)
MFEM_DEPRECATED HypreParMatrix * GetData()
Deprecated. Use HypreTriSolve::GetData() const instead.
void SetMassMatrix(const HypreParMatrix &M)
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
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.
void HypreRead() const
Prepare the HypreParVector for read access in hypre's device memory space, HYPRE_MEMORY_DEVICE.
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.
void HypreRead() const
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
HYPRE_BigInt GlobalSize() const
Returns the global number of rows.
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.
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.
HypreParMatrix & Add(const double beta, const HypreParMatrix &B)
HypreParVector * B
Right-hand side and solution vector.
Wrapper for Hypre's native parallel ILU preconditioner.
void GetNumIterations(int &num_iterations) const
const HypreParMatrix * GetData() const
double window_params[3]
Parameters for windowing function of FIR filter.
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
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.
void SetPreconditioner(HypreSolver &precond)
HypreILU()
Constructor; sets the default options.
void SetRestriction(int restrict_type)
Expert option - consult hypre documentation/team.
HypreAMS(ParFiniteElementSpace *edge_fespace)
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 SetMassMatrix(Operator &M)
Abstract parallel finite element space.
HypreParVector()
Default constructor, no underlying hypre_ParVector is created.
void SetPrintLevel(int logging)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void SetAdvectiveOptions(int distance=15, const std::string &prerelax="", const std::string &postrelax="FFC")
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
GMRES Solve 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.
HypreParMatrix & operator+=(const HypreParMatrix &B)
bool pos_l1_norms
If set, take absolute values of the computed l1_norms.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
void SetPrintLevel(int print_lvl)
void SetZeroInitialIterate()
non-hypre setting
void SetStrongThresholdR(double strengthR)
Expert option - consult hypre documentation/team.
MPI_Comm GetComm() const
MPI communicator.
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
virtual MFEM_DEPRECATED N_Vector ToNVector()
(DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_PARALLEL.
double poly_fraction
Fraction of spectrum to smooth for polynomial relaxation.
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...
void GetEigenvalues(Array< double > &eigenvalues) const
Collect the converged eigenvalues.
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...
HypreFGMRES(MPI_Comm comm)
HYPRE_BigInt N() const
Returns the global number of columns.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's FGMRES.
HYPRE_BigInt M() const
Returns the global number of rows.
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.
HypreGMRES(MPI_Comm comm)
void GetOffd(SparseMatrix &offd, HYPRE_BigInt *&cmap) const
Get the local off-diagonal block. NOTE: 'offd' will not own any data.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's GMRES.
virtual void AssembleDiagonal(Vector &diag) const
Return the diagonal of the matrix (Operator interface).
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
void Print(const char *fname) const
Prints the locally owned rows in parallel.
void HypreWrite()
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.
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...
HypreLOBPCG(MPI_Comm comm)
HYPRE_BigInt * GetRowStarts() const
Return the parallel row partitioning array.
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...
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=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 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 ...
virtual ~HypreParaSails()
The ParaSails preconditioner in hypre.
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.
void PrintHash(std::ostream &out) const
Print sizes and hashes for all data arrays of the HypreParMatrix from the local MPI rank...
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)
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)
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
FGMRES Setup function.
void SetOperatorSymmetry(bool is_sym)
void GetBlocks(Array2D< HypreParMatrix * > &blocks, bool interleaved_rows=false, bool interleaved_cols=false) const
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
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 SetTol(double tol)
Expert option - consult hypre documentation/team.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
PCG Solve function.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
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.
HypreParMatrix * Transpose() const
Returns the transpose of *this.
HYPRE_Int Randomize(HYPRE_Int seed)
Set random values.
void SetFilterThresholdR(double filterR)
Expert option - consult hypre documentation/team.
Ignore hypre errors (see e.g. HypreADS)
double relax_weight
Damping coefficient (usually <= 1)
Parallel smoothers in hypre.
void SetElasticityOptions(ParFiniteElementSpace *fespace)
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
ILU Solve function.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
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.
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
hypre's internal Solve function
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)
HYPRE_BigInt GetGlobalNumCols() const
Return the global number of columns.
HypreParMatrix()
An empty matrix to be used as a reference to an existing matrix.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
bool A_is_symmetric
A flag that indicates whether the linear system matrix A is symmetric.
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)
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
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)
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
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...
signed char OwnsOffd() const
Get offd ownership flag.
void HostRead() const
Update the internal hypre_ParCSRMatrix object, A, to be on host.
MemoryType
Memory types supported by MFEM.
signed char OwnsColMap() const
Get colmap ownership flag.
void SetPrintLevel(int print_lvl)
void Solve()
Solve the eigenproblem.
void SetData(double *data_)
Sets the data of the Vector and the hypre_ParVector to data_.
void SetNumModes(int num_eigs)
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
HypreParMatrix(hypre_ParCSRMatrix *a, bool owner=true)
Converts hypre's format to HypreParMatrix.
virtual ~HypreDiagScale()
const HypreParVector & GetEigenvector(unsigned int i) const
Extract a single eigenvector.
int GetNumRows() const
Returns the number of rows in the diagonal block of the ParCSRMatrix.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
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.
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)
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 GetNumIterations(int &num_iterations) const
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
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.
signed char OwnsDiag() const
Get diag ownership flag.
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 Print(const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0) const
Prints the locally owned rows in parallel.
void SetPrintLevel(HYPRE_Int print_level)
Set the print level: 0 = none, 1 = setup, 2 = solve, 3 = setup+solve.
Host memory; using new[] and delete[].
void SetSubSpaceProjector(Operator &proj)
virtual MFEM_DEPRECATED N_Vector ToNVector()
(DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_SERIAL. ...
constexpr MemoryClass GetHypreMemoryClass()
The MemoryClass used by Hypre objects.
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
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.
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.
virtual void MultTranspose(const Vector &b, Vector &x) const
Apply transpose of the smoother to relax the linear system Ax=b.
int GetNumCols() const
Returns the number of columns in the diagonal block of the ParCSRMatrix.
void SetZeroInitialIterate()
non-hypre setting
void SetOperator(Operator &A)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
int relax_times
Number of relaxation sweeps.
void GetEigenvalues(Array< double > &eigenvalues) const
Collect the converged eigenvalues.
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.
HYPRE_BigInt NNZ() const
Returns the global number of nonzeros.
const HypreParVector & GetEigenvector(unsigned int i) const
Extract a single eigenvector.
HYPRE_BigInt * GetColStarts() const
Return the parallel column partitioning array.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
ILU Setup function.
const HYPRE_BigInt * Partitioning() const
Returns the parallel row/column partitioning.
HypreParVector & operator=(double d)
Set constant values.
void SetAbsTol(double tol)
void PrintCommPkg(std::ostream &out=mfem::out) const
Print information about the hypre_ParCSRCommPkg of the HypreParMatrix.
HypreParMatrix & operator=(double value)
Initialize all entries with value.
HypreParMatrix * ExtractSubmatrix(const Array< int > &indices, double threshold=0.0) const
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.
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A * x + beta * y.
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)
const HYPRE_BigInt * RowPart() const
Returns the row partitioning (const version)
double omega
SOR parameter (usually in (0,2))
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 HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
void SetNodal(int blocksize)
Expert option - consult hypre documentation/team.
int GetOwnership() const
Gets ownership of the internal hypre_ParVector.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
PCG Setup function.
void WrapHypreParCSRMatrix(hypre_ParCSRMatrix *a, bool owner=true)
Converts hypre's format to HypreParMatrix.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Wrapper for hypre's ParCSR matrix class.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
void HostReadWrite()
Update the internal hypre_ParCSRMatrix object, A, to be on host.
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.
Memory< double > auxB
Auxiliary buffers for the case when the input or output arrays in methods like Mult(const Vector &...
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
~HypreParVector()
Calls hypre's destroy function.
void SetNumModes(int num_eigs)
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
HypreParMatrix * HypreParMatrixFromBlocks(Array2D< HypreParMatrix * > &blocks, Array2D< double > *blockCoeff)
Returns a merged hypre matrix constructed from hypre matrix blocks.
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.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
MemoryClass
Memory classes identify sets of memory types.
void WrapMemoryReadWrite(Memory< double > &mem)
Replace the HypreParVector's data with the given Memory, mem, and prepare the vector for read and wri...
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
GMRES Setup function.
int poly_order
Order of the smoothing polynomial.
double lambda
Taubin's lambda-mu method parameters.
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()
MFEM_DEPRECATED HypreParMatrix * GetData()
Deprecated. Use HypreDiagScale::GetData() const instead.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
void HypreReadWrite()
Update the internal hypre_ParCSRMatrix object, A, to be in hypre memory space.