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"
37 #ifdef MFEM_USE_SUNDIALS
38 #include <nvector/nvector_parhyp.h>
44 class ParFiniteElementSpace;
50 template <
typename int_type>
51 inline int to_int(int_type i)
53 MFEM_ASSERT(int_type(
int(i)) == i,
"overflow converting int_type to int");
58 template <>
inline int to_int(
int i) {
return i; }
63 inline int to_int(HYPRE_Int i)
65 MFEM_ASSERT(HYPRE_Int(
int(i)) == i,
"overflow converting HYPRE_Int to int");
84 inline void _SetDataAndSize_();
89 HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col);
115 operator hypre_ParVector*()
const {
return x; }
116 #ifndef HYPRE_PAR_VECTOR_STRUCT
118 operator HYPRE_ParVector()
const {
return (HYPRE_ParVector) x; }
148 void Print(
const char *fname)
const;
153 #ifdef MFEM_USE_SUNDIALS
156 virtual N_Vector
ToNVector() {
return N_VMake_ParHyp(x); }
165 double InnerProduct(HypreParVector &x, HypreParVector &y);
166 double InnerProduct(HypreParVector *x, HypreParVector *y);
171 double ParNormlp(
const Vector &vec,
double p, MPI_Comm comm);
179 hypre_ParCSRMatrix *A;
198 signed char diagOwner, offdOwner, colMapOwner;
201 signed char ParCSROwner;
212 static char CopyCSR(
SparseMatrix *csr, hypre_CSRMatrix *hypre_csr);
216 static char CopyBoolCSR(
Table *bool_csr, hypre_CSRMatrix *hypre_csr);
220 static void CopyCSR_J(hypre_CSRMatrix *hypre_csr,
int *J);
232 if (!owner) { ParCSROwner = 0; }
243 HypreParMatrix(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *row_starts,
250 HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
256 HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
263 HYPRE_Int global_num_rows, HYPRE_Int global_num_cols,
264 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
265 HYPRE_Int *diag_i, HYPRE_Int *diag_j,
double *diag_data,
266 HYPRE_Int *offd_i, HYPRE_Int *offd_j,
double *offd_data,
267 HYPRE_Int offd_num_cols, HYPRE_Int *offd_col_map);
270 HypreParMatrix(MPI_Comm comm, HYPRE_Int *row_starts, HYPRE_Int *col_starts,
276 HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
277 HYPRE_Int *col_starts,
Table *diag);
282 HypreParMatrix(MPI_Comm comm,
int id,
int np, HYPRE_Int *row, HYPRE_Int *col,
283 HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd,
284 HYPRE_Int *j_offd, HYPRE_Int *cmap, HYPRE_Int cmap_size);
291 HYPRE_Int glob_ncols,
int *I, HYPRE_Int *J,
292 double *data, HYPRE_Int *rows, HYPRE_Int *cols);
301 operator hypre_ParCSRMatrix*()
const {
return A; }
302 #ifndef HYPRE_PAR_CSR_MATRIX_STRUCT
304 operator HYPRE_ParCSRMatrix() {
return (HYPRE_ParCSRMatrix) A; }
311 { diagOwner = diag, offdOwner = offd, colMapOwner = colmap; }
330 inline HYPRE_Int
NNZ()
const {
return A->num_nonzeros; }
332 inline HYPRE_Int *
RowPart() {
return A->row_starts; }
334 inline HYPRE_Int *
ColPart() {
return A->col_starts; }
336 inline const HYPRE_Int *
RowPart()
const {
return A->row_starts; }
338 inline const HYPRE_Int *
ColPart()
const {
return A->col_starts; }
340 inline HYPRE_Int
M()
const {
return A->global_num_rows; }
342 inline HYPRE_Int
N()
const {
return A->global_num_cols; }
354 bool interleaved_rows =
false,
355 bool interleaved_cols =
false)
const;
364 hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(A)));
371 hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(A)));
375 {
return hypre_ParCSRMatrixGlobalNumRows(A); }
378 {
return hypre_ParCSRMatrixGlobalNumCols(A); }
380 HYPRE_Int *
GetRowStarts()
const {
return hypre_ParCSRMatrixRowStarts(A); }
382 HYPRE_Int *
GetColStarts()
const {
return hypre_ParCSRMatrixColStarts(A); }
386 double alpha = 1.0,
double beta = 0.0);
388 HYPRE_Int
Mult(HYPRE_ParVector x, HYPRE_ParVector y,
389 double alpha = 1.0,
double beta = 0.0);
392 double alpha = 1.0,
double beta = 0.0);
398 {
Mult(1.0, x, 0.0, y); }
406 internal::hypre_ParCSRMatrixBooleanMatvec(A, alpha, x, beta, y);
411 { internal::hypre_ParCSRMatrixSetConstantValues(A, value);
return *
this; }
423 MFEM_VERIFY(internal::hypre_ParCSRMatrixSum(A, beta, B.A) == 0,
424 "error in hypre_ParCSRMatrixSum");
436 HYPRE_Int* row_starts = NULL)
const;
462 void Print(
const char *fname, HYPRE_Int offi = 0, HYPRE_Int offj = 0);
464 void Read(MPI_Comm comm,
const char *fname);
480 HypreParMatrix *
Add(
double alpha,
const HypreParMatrix &A,
481 double beta,
const HypreParMatrix &B);
484 HypreParMatrix *
ParMult(
const HypreParMatrix *A,
const HypreParMatrix *B);
488 HypreParMatrix *
ParAdd(
const HypreParMatrix *A,
const HypreParMatrix *B);
491 HypreParMatrix *
RAP(
const HypreParMatrix *A,
const HypreParMatrix *P);
493 HypreParMatrix *
RAP(
const HypreParMatrix * Rt,
const HypreParMatrix *A,
494 const HypreParMatrix *P);
499 void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae,
500 const Array<int> &ess_dof_list,
const Vector &X, Vector &B);
625 virtual operator HYPRE_Solver()
const = 0;
628 virtual HYPRE_PtrToParSolverFcn
SetupFcn()
const = 0;
630 virtual HYPRE_PtrToParSolverFcn
SolveFcn()
const = 0;
633 {
mfem_error(
"HypreSolvers do not support SetOperator!"); }
646 HYPRE_Solver pcg_solver;
670 HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_it);
675 virtual operator HYPRE_Solver()
const {
return pcg_solver; }
679 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSetup; }
682 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSolve; }
695 HYPRE_Solver gmres_solver;
713 virtual operator HYPRE_Solver()
const {
return gmres_solver; }
717 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSetup; }
720 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSolve; }
733 virtual operator HYPRE_Solver()
const {
return NULL; }
736 {
return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentitySetup; }
738 {
return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentity; }
749 virtual operator HYPRE_Solver()
const {
return NULL; }
752 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScaleSetup; }
754 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScale; }
764 HYPRE_Solver sai_precond;
772 virtual operator HYPRE_Solver()
const {
return sai_precond; }
775 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSetup; }
777 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSolve; }
786 HYPRE_Solver amg_precond;
795 void RecomputeRBMs();
798 void SetDefaultOptions();
803 void ResetAMGPrecond();
825 { HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level); }
828 virtual operator HYPRE_Solver()
const {
return amg_precond; }
831 {
return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSetup; }
833 {
return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSolve; }
839 HypreParMatrix*
DiscreteGrad(ParFiniteElementSpace *edge_fespace,
840 ParFiniteElementSpace *vert_fespace);
842 HypreParMatrix*
DiscreteCurl(ParFiniteElementSpace *face_fespace,
843 ParFiniteElementSpace *edge_fespace);
867 virtual operator HYPRE_Solver()
const {
return ams; }
870 {
return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSetup; }
872 {
return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSolve; }
900 virtual operator HYPRE_Solver()
const {
return ads; }
903 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSetup; }
905 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSolve; }
945 HYPRE_Solver lobpcg_solver;
948 mv_InterfaceInterpreter interpreter;
951 HYPRE_MatvecFunctions matvec_fn;
957 class HypreMultiVector;
960 HypreMultiVector * multi_vec;
969 class HypreMultiVector
973 mv_MultiVectorPtr mv_ptr;
983 mv_InterfaceInterpreter & interpreter);
987 void Randomize(HYPRE_Int seed);
995 operator mv_MultiVectorPtr()
const {
return mv_ptr; }
997 mv_MultiVectorPtr & GetMultiVector() {
return mv_ptr; }
1000 static void * OperatorMatvecCreate(
void *A,
void *x );
1001 static HYPRE_Int OperatorMatvec(
void *matvec_data,
1002 HYPRE_Complex
alpha,
1007 static HYPRE_Int OperatorMatvecDestroy(
void *matvec_data );
1009 static HYPRE_Int PrecondSolve(
void *solver,
1013 static HYPRE_Int PrecondSetup(
void *solver,
1081 HYPRE_Solver ame_solver;
1087 HYPRE_Real * eigenvalues;
1090 HYPRE_ParVector * multi_vec;
1094 void createDummyVectors();
1126 #endif // MFEM_USE_MPI
virtual ~HypreBoomerAMG()
const HYPRE_Int * ColPart() const
Returns the column partitioning (const version)
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)
MPI_Comm GetComm() const
MPI communicator.
Vector * GlobalVector() const
Returns the global vector in each processor.
The Auxiliary-space Maxwell Solver in hypre.
HypreParVector * X0
FIR Filter Temporary Vectors.
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
double min_eig_est
Minimal eigenvalue estimate for polynomial smoothing.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
The Auxiliary-space Divergence Solver in hypre.
void SetPrintLevel(int print_lvl)
void Print(const char *fname, HYPRE_Int offi=0, HYPRE_Int offj=0)
Prints the locally owned rows in parallel.
MPI_Comm GetComm()
MPI communicator.
HYPRE_Int N() const
Returns the global number of columns.
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
int setup_called
Was hypre's Setup function called already?
HYPRE_Int * ColPart()
Returns the column partitioning.
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 * LeftDiagMult(const SparseMatrix &D, HYPRE_Int *row_starts=NULL) const
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A^t * x + beta * y.
HypreParMatrix & Add(const double beta, const HypreParMatrix &B)
HypreParVector * B
Right-hand side and solution vector.
HypreDiagScale(HypreParMatrix &A)
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)
HypreADS(HypreParMatrix &A, ParFiniteElementSpace *face_fespace)
void SetWindowByName(const char *window_name)
Convenience function for setting canonical windowing parameters.
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
void SetPreconditioner(HypreSolver &precond)
virtual N_Vector ToNVector()
Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_PARHYP.
void SetPreconditioner(Solver &precond)
void SetMassMatrix(Operator &M)
Abstract parallel finite element space.
void SetPrintLevel(int logging)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void SetOperator(HypreParMatrix &A)
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
GMRES Solve function.
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
HypreParMatrix & operator+=(const HypreParMatrix &B)
bool pos_l1_norms
If set, take absolute values of the computed l1_norms.
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
void SetPrintLevel(int print_lvl)
void ScaleRows(const Vector &s)
Scale the local row i by s(i).
double poly_fraction
Fraction of spectrum to smooth for polynomial relaxation.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A * B.
void SetPrintLevel(int print_lvl)
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}.
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.
HYPRE_Int NNZ() const
Returns the global number of nonzeros.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's GMRES.
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 EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae, const Array< int > &ess_dof_list, const Vector &X, Vector &B)
void SetSystemsOptions(int dim)
HypreLOBPCG(MPI_Comm comm)
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Relax the linear system Ax=b.
hypre_ParVector * StealParVector()
Changes the ownership of the the vector.
The BoomerAMG solver in hypre.
HYPRE_Int GetGlobalNumRows() const
void SetSymmetry(int sym)
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.
HYPRE_Int GetGlobalNumCols() const
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 SetLogging(int logging)
void SetZeroInintialIterate()
non-hypre setting
Jacobi preconditioner in hypre.
HyprePCG(HypreParMatrix &_A)
virtual HYPRE_PtrToParSolverFcn SolveFcn() const =0
hypre's internal Solve function
void SetRelTol(double rel_tol)
void SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0)
int to_int(const std::string &str)
void SetPrintLevel(int print_level)
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.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
PCG Solve function.
HypreParVector ** StealEigenvectors()
Transfer ownership of the converged eigenvectors.
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.
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 SetupFcn() const
hypre's internal Setup function
void BooleanMult(int alpha, int *x, int beta, int *y)
HypreParMatrix * A
The linear system matrix.
HypreAMS(HypreParMatrix &A, ParFiniteElementSpace *edge_fespace)
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
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.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
void SetMaxIter(int max_iter)
HypreParMatrix * GetData()
HYPRE_Int * RowPart()
Returns the row partitioning.
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
Wrapper for hypre's parallel vector class.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
Type
Enumeration defining IDs for some classes derived from Operator.
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
signed char OwnsOffd() const
Get offd ownership flag.
const HYPRE_Int * RowPart() const
Returns the row partitioning (const version)
signed char OwnsColMap() const
Get colmap ownership flag.
void SetPrintLevel(int print_lvl)
void mfem_error(const char *msg)
void Solve()
Solve the eigenproblem.
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()
int GetNumRows() const
Returns the number of rows in the diagonal block of the ParCSRMatrix.
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.
void SetRelTol(double rel_tol)
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
int height
Dimension of the output / number of rows in the matrix.
void SetMaxIter(int max_iter)
void SetFIRCoefficients(double max_eig)
Compute window and Chebyshev coefficients for given polynomial order.
void SetRandomSeed(int s)
signed char OwnsDiag() const
Get diag ownership flag.
double * fir_coeffs
Combined coefficients for windowing and Chebyshev polynomials.
void Threshold(double threshold=0.0)
Remove values smaller in absolute value than some threshold.
double InnerProduct(HypreParVector *x, HypreParVector *y)
HYPRE_Int GlobalSize()
Returns the global number of rows.
HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col)
void SetSubSpaceProjector(Operator &proj)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void SetData(double *_data)
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
void SetTaubinOptions(double lambda, double mu, int iter)
Set parameters for Taubin's lambda-mu method.
HypreParVector * B
Right-hand side and solution vectors.
void SetMassMatrix(HypreParMatrix &M)
HYPRE_Int M() const
Returns the global number of rows.
void SetZeroInintialIterate()
non-hypre setting
int GetNumCols() const
Returns the number of columns in the diagonal block of the ParCSRMatrix.
void SetOperator(Operator &A)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
void SetPolyOptions(int poly_order, double poly_fraction)
Set parameters for polynomial smoothing.
void GetNumIterations(int &num_iterations)
int relax_times
Number of relaxation sweeps.
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...
void SetOwnership(int own)
Sets ownership of the internal hypre_ParVector.
Abstract class for hypre's solvers and preconditioners.
HypreParVector & operator=(double d)
Set constant values.
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.
HypreGMRES(HypreParMatrix &_A)
hypre_ParCSRMatrix * StealData()
Changes the ownership of the the matrix.
void SetPrecondUsageMode(int pcg_mode)
virtual ~HypreParMatrix()
Calls hypre's destroy function.
HypreParMatrix * A
The linear system matrix.
void SetMaxIter(int max_iter)
void GetOffd(SparseMatrix &offd, HYPRE_Int *&cmap) const
Get the local off-diagonal block. NOTE: 'offd' will not own any data.
HYPRE_Int * GetColStarts() const
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
double omega
SOR parameter (usually in (0,2))
HYPRE_Int * GetRowStarts() const
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
int GetOwnership() const
Gets ownership of the internal hypre_ParVector.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
PCG Setup function.
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 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.
void SetSingularProblem()
Set this option when solving a curl-curl problem with zero mass term.
void Solve()
Solve the eigenproblem.
void Read(MPI_Comm comm, const char *fname)
Reads the matrix from a file.
HypreParVector * V
Temporary vectors.
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
HypreParaSails(HypreParMatrix &A)
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
GMRES Setup function.
int width
Dimension of the input / number of columns in the matrix.
int poly_order
Order of the smoothing polynomial.
HYPRE_Int * Partitioning()
Returns the row partitioning.
double lambda
Taubin's lambda-mu method parameters.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function