15 #include "../config/config.hpp"
26 #include "_hypre_parcsr_mv.h"
27 #include "_hypre_parcsr_ls.h"
28 #include "temp_multivector.h"
31 #error "MFEM does not work with HYPRE's complex numbers support"
40 class ParFiniteElementSpace;
47 inline int to_int(HYPRE_Int i)
50 MFEM_ASSERT(HYPRE_Int(
int(i)) == i,
"overflow converting HYPRE_Int to int");
69 inline void _SetDataAndSize_();
74 HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col);
98 operator hypre_ParVector*()
const;
99 #ifndef HYPRE_PAR_VECTOR_STRUCT
100 operator HYPRE_ParVector()
const;
103 hypre_ParVector *
StealParVector() { own_ParVector = 0;
return x; }
131 void Print(
const char *fname)
const;
138 double InnerProduct(HypreParVector &x, HypreParVector &y);
139 double InnerProduct(HypreParVector *x, HypreParVector *y);
147 hypre_ParCSRMatrix *A;
166 char diagOwner, offdOwner, colMapOwner;
180 static char CopyCSR(
SparseMatrix *csr, hypre_CSRMatrix *hypre_csr);
184 static char CopyBoolCSR(
Table *bool_csr, hypre_CSRMatrix *hypre_csr);
188 static void CopyCSR_J(hypre_CSRMatrix *hypre_csr,
int *J);
206 HypreParMatrix(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *row_starts,
213 HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
219 HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
226 HYPRE_Int global_num_rows, HYPRE_Int global_num_cols,
227 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
228 HYPRE_Int *diag_i, HYPRE_Int *diag_j,
double *diag_data,
229 HYPRE_Int *offd_i, HYPRE_Int *offd_j,
double *offd_data,
230 HYPRE_Int offd_num_cols, HYPRE_Int *offd_col_map);
233 HypreParMatrix(MPI_Comm comm, HYPRE_Int *row_starts, HYPRE_Int *col_starts,
239 HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
240 HYPRE_Int *col_starts,
Table *diag);
245 HypreParMatrix(MPI_Comm comm,
int id,
int np, HYPRE_Int *row, HYPRE_Int *col,
246 HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd,
247 HYPRE_Int *j_offd, HYPRE_Int *cmap, HYPRE_Int cmap_size);
254 HYPRE_Int glob_ncols,
int *I, HYPRE_Int *J,
255 double *data, HYPRE_Int *rows, HYPRE_Int *cols);
264 operator hypre_ParCSRMatrix*() {
return A; }
265 #ifndef HYPRE_PAR_CSR_MATRIX_STRUCT
266 operator HYPRE_ParCSRMatrix() {
return (HYPRE_ParCSRMatrix) A; }
274 { diagOwner = diag, offdOwner = offd, colMapOwner = colmap; }
293 inline HYPRE_Int
NNZ() {
return A->num_nonzeros; }
295 inline HYPRE_Int *
RowPart() {
return A->row_starts; }
297 inline HYPRE_Int *
ColPart() {
return A->col_starts; }
299 inline HYPRE_Int
M() {
return A->global_num_rows; }
301 inline HYPRE_Int
N() {
return A->global_num_cols; }
313 bool interleaved_rows =
false,
314 bool interleaved_cols =
false)
const;
322 return internal::to_int(
323 hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(A)));
329 return internal::to_int(
330 hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(A)));
334 {
return hypre_ParCSRMatrixGlobalNumRows(A); }
337 {
return hypre_ParCSRMatrixGlobalNumCols(A); }
339 HYPRE_Int *
GetRowStarts()
const {
return hypre_ParCSRMatrixRowStarts(A); }
341 HYPRE_Int *
GetColStarts()
const {
return hypre_ParCSRMatrixColStarts(A); }
345 double alpha = 1.0,
double beta = 0.0);
347 HYPRE_Int
Mult(HYPRE_ParVector x, HYPRE_ParVector y,
348 double alpha = 1.0,
double beta = 0.0);
351 double alpha = 1.0,
double beta = 0.0);
357 {
Mult(1.0, x, 0.0, y); }
365 internal::hypre_ParCSRMatrixBooleanMatvec(A, alpha, x, beta, y);
376 HYPRE_Int* row_starts = NULL)
const;
402 void Print(
const char *fname, HYPRE_Int offi = 0, HYPRE_Int offj = 0);
404 void Read(MPI_Comm comm,
const char *fname);
411 HypreParMatrix *
ParMult(HypreParMatrix *A, HypreParMatrix *B);
414 HypreParMatrix *
RAP(HypreParMatrix *A, HypreParMatrix *P);
416 HypreParMatrix *
RAP(HypreParMatrix * Rt, HypreParMatrix *A, HypreParMatrix *P);
421 void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae,
422 const Array<int> &ess_dof_list,
const Vector &X, Vector &B);
540 virtual operator HYPRE_Solver()
const = 0;
543 virtual HYPRE_PtrToParSolverFcn
SetupFcn()
const = 0;
545 virtual HYPRE_PtrToParSolverFcn
SolveFcn()
const = 0;
548 {
mfem_error(
"HypreSolvers do not support SetOperator!"); }
561 HYPRE_Solver pcg_solver;
585 HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_it);
586 num_iterations = internal::to_int(num_it);
590 virtual operator HYPRE_Solver()
const {
return pcg_solver; }
594 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSetup; }
597 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSolve; }
610 HYPRE_Solver gmres_solver;
628 virtual operator HYPRE_Solver()
const {
return gmres_solver; }
632 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSetup; }
635 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSolve; }
648 virtual operator HYPRE_Solver()
const {
return NULL; }
651 {
return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentitySetup; }
653 {
return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentity; }
664 virtual operator HYPRE_Solver()
const {
return NULL; }
667 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScaleSetup; }
669 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScale; }
679 HYPRE_Solver sai_precond;
687 virtual operator HYPRE_Solver()
const {
return sai_precond; }
690 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSetup; }
692 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSolve; }
701 HYPRE_Solver amg_precond;
710 void RecomputeRBMs();
713 void SetDefaultOptions();
718 void ResetAMGPrecond();
740 { HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level); }
743 virtual operator HYPRE_Solver()
const {
return amg_precond; }
746 {
return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSetup; }
748 {
return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSolve; }
754 HypreParMatrix*
DiscreteGrad(ParFiniteElementSpace *edge_fespace,
755 ParFiniteElementSpace *vert_fespace);
757 HypreParMatrix*
DiscreteCurl(ParFiniteElementSpace *face_fespace,
758 ParFiniteElementSpace *edge_fespace);
782 virtual operator HYPRE_Solver()
const {
return ams; }
785 {
return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSetup; }
787 {
return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSolve; }
815 virtual operator HYPRE_Solver()
const {
return ads; }
818 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSetup; }
820 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSolve; }
860 HYPRE_Solver lobpcg_solver;
863 mv_InterfaceInterpreter interpreter;
866 HYPRE_MatvecFunctions matvec_fn;
872 class HypreMultiVector;
875 HypreMultiVector * multi_vec;
881 class HypreMultiVector
885 mv_MultiVectorPtr mv_ptr;
895 mv_InterfaceInterpreter & interpreter);
899 void Randomize(HYPRE_Int seed);
907 operator mv_MultiVectorPtr()
const {
return mv_ptr; }
909 mv_MultiVectorPtr & GetMultiVector() {
return mv_ptr; }
912 static void * OperatorMatvecCreate(
void *A,
void *x );
913 static HYPRE_Int OperatorMatvec(
void *matvec_data,
919 static HYPRE_Int OperatorMatvecDestroy(
void *matvec_data );
921 static HYPRE_Int PrecondSolve(
void *solver,
925 static HYPRE_Int PrecondSetup(
void *solver,
990 HYPRE_Solver ame_solver;
996 HYPRE_Real * eigenvalues;
999 HYPRE_ParVector * multi_vec;
1003 void createDummyVectors();
1034 #endif // MFEM_USE_MPI
virtual ~HypreBoomerAMG()
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
HYPRE_Int NNZ()
Returns the global number of nonzeros.
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
void SetOwnerFlags(char diag, char offd, char colmap)
Explicitly set the three ownership flags, see docs for diagOwner etc.
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 M()
Returns the global number of rows.
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'.
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.
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.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
HypreADS(HypreParMatrix &A, ParFiniteElementSpace *face_fespace)
void SetWindowByName(const char *window_name)
Convenience function for setting canonical windowing parameters.
HypreParMatrix(hypre_ParCSRMatrix *a)
Converts hypre's format to HypreParMatrix.
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
void SetPreconditioner(HypreSolver &precond)
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.
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.
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.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's GMRES.
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)
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.
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 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 SetResidualConvergenceOptions(int res_frequency=-1, double rtol=0.0)
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.
HypreParMatrix * RAP(HypreParMatrix *A, HypreParMatrix *P)
Returns the matrix P^t * A * P.
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.
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)
char OwnsColMap() const
Get colmap ownership flag.
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
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
char OwnsOffd() const
Get offd 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
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.
double max_eig_est
Maximal eigenvalue estimate for polynomial smoothing.
HypreParMatrix * Transpose()
Returns the transpose of *this.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
hypre's internal Setup function
void SetMaxIter(int max_iter)
void SetFIRCoefficients(double max_eig)
Compute window and Chebyshev coefficients for given polynomial order.
void SetRandomSeed(int s)
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)
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)
HYPRE_Int N()
Returns the global number of columns.
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)
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)
void SetOwnership(int own)
Sets ownership of the internal hypre_ParVector.
Abstract class for hypre's solvers and preconditioners.
char OwnsDiag() const
Get diag ownership flag.
HypreParVector & operator=(double d)
Set constant values.
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.
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 poly_order
Order of the smoothing polynomial.
HYPRE_Int * Partitioning()
Returns the row partitioning.
double lambda
Taubin's lambda-mu method parameters.
HypreParMatrix * ParMult(HypreParMatrix *A, HypreParMatrix *B)
Returns the matrix A * B.
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function