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"
41 class ParFiniteElementSpace;
47 template <
typename int_type>
48 inline int to_int(int_type i)
50 MFEM_ASSERT(int_type(
int(i)) == i,
"overflow converting int_type to int");
55 template <>
inline int to_int(
int i) {
return i; }
60 inline int to_int(HYPRE_Int i)
62 MFEM_ASSERT(HYPRE_Int(
int(i)) == i,
"overflow converting HYPRE_Int to int");
81 inline void _SetDataAndSize_();
96 HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col);
125 operator hypre_ParVector*()
const {
return x; }
126 #ifndef HYPRE_PAR_VECTOR_STRUCT
128 operator HYPRE_ParVector()
const {
return (HYPRE_ParVector) x; }
157 void Print(
const char *fname)
const;
162 #ifdef MFEM_USE_SUNDIALS
165 MFEM_DEPRECATED
virtual N_Vector
ToNVector();
171 double InnerProduct(HypreParVector &x, HypreParVector &y);
172 double InnerProduct(HypreParVector *x, HypreParVector *y);
177 double ParNormlp(
const Vector &vec,
double p, MPI_Comm comm);
185 hypre_ParCSRMatrix *A;
204 signed char diagOwner, offdOwner, colMapOwner;
207 signed char ParCSROwner;
217 static char CopyCSR(
SparseMatrix *csr, hypre_CSRMatrix *hypre_csr);
221 static char CopyBoolCSR(
Table *bool_csr, hypre_CSRMatrix *hypre_csr);
225 static void CopyCSR_J(hypre_CSRMatrix *hypre_csr,
int *J);
237 if (!owner) { ParCSROwner = 0; }
251 HypreParMatrix(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *row_starts,
260 HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
261 HYPRE_Int *col_starts,
269 HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
279 HYPRE_Int global_num_rows, HYPRE_Int global_num_cols,
280 HYPRE_Int *row_starts, HYPRE_Int *col_starts,
281 HYPRE_Int *diag_i, HYPRE_Int *diag_j,
double *diag_data,
282 HYPRE_Int *offd_i, HYPRE_Int *offd_j,
double *offd_data,
283 HYPRE_Int offd_num_cols,
284 HYPRE_Int *offd_col_map);
289 HypreParMatrix(MPI_Comm comm, HYPRE_Int *row_starts, HYPRE_Int *col_starts,
297 HYPRE_Int global_num_cols, HYPRE_Int *row_starts,
298 HYPRE_Int *col_starts,
306 HypreParMatrix(MPI_Comm comm,
int id,
int np, HYPRE_Int *row, HYPRE_Int *col,
307 HYPRE_Int *i_diag, HYPRE_Int *j_diag, HYPRE_Int *i_offd,
308 HYPRE_Int *j_offd, HYPRE_Int *cmap,
309 HYPRE_Int cmap_size);
318 HYPRE_Int glob_ncols,
int *I, HYPRE_Int *J,
319 double *data, HYPRE_Int *rows,
333 operator hypre_ParCSRMatrix*()
const {
return A; }
334 #ifndef HYPRE_PAR_CSR_MATRIX_STRUCT
336 operator HYPRE_ParCSRMatrix() {
return (HYPRE_ParCSRMatrix) A; }
343 { diagOwner = diag, offdOwner = offd, colMapOwner = colmap; }
362 inline HYPRE_Int
NNZ()
const {
return A->num_nonzeros; }
366 inline HYPRE_Int *
RowPart() {
return A->row_starts; }
370 inline HYPRE_Int *
ColPart() {
return A->col_starts; }
374 inline const HYPRE_Int *
RowPart()
const {
return A->row_starts; }
378 inline const HYPRE_Int *
ColPart()
const {
return A->col_starts; }
380 inline HYPRE_Int
M()
const {
return A->global_num_rows; }
382 inline HYPRE_Int
N()
const {
return A->global_num_cols; }
394 bool interleaved_rows =
false,
395 bool interleaved_cols =
false)
const;
404 hypre_CSRMatrixNumRows(hypre_ParCSRMatrixDiag(A)));
411 hypre_CSRMatrixNumCols(hypre_ParCSRMatrixDiag(A)));
416 {
return hypre_ParCSRMatrixGlobalNumRows(A); }
420 {
return hypre_ParCSRMatrixGlobalNumCols(A); }
425 HYPRE_Int *
GetRowStarts()
const {
return hypre_ParCSRMatrixRowStarts(A); }
430 HYPRE_Int *
GetColStarts()
const {
return hypre_ParCSRMatrixColStarts(A); }
434 double alpha = 1.0,
double beta = 0.0);
436 HYPRE_Int
Mult(HYPRE_ParVector x, HYPRE_ParVector y,
437 double alpha = 1.0,
double beta = 0.0);
440 double alpha = 1.0,
double beta = 0.0);
446 {
Mult(1.0, x, 0.0, y); }
460 internal::hypre_ParCSRMatrixBooleanMatvec(A, alpha, const_cast<int*>(x),
468 internal::hypre_ParCSRMatrixBooleanMatvecT(A, alpha, const_cast<int*>(x),
474 { internal::hypre_ParCSRMatrixSetConstantValues(A, value);
return *
this; }
489 MFEM_VERIFY(internal::hypre_ParCSRMatrixSum(A, beta, B.A) == 0,
490 "error in hypre_ParCSRMatrixSum");
503 HYPRE_Int* row_starts = NULL)
const;
537 void Print(
const char *fname, HYPRE_Int offi = 0, HYPRE_Int offj = 0);
539 void Read(MPI_Comm comm,
const char *fname);
555 HypreParMatrix *
Add(
double alpha,
const HypreParMatrix &A,
556 double beta,
const HypreParMatrix &B);
560 HypreParMatrix *
ParMult(
const HypreParMatrix *A,
const HypreParMatrix *B,
561 bool own_matrix =
false);
565 HypreParMatrix *
ParAdd(
const HypreParMatrix *A,
const HypreParMatrix *B);
568 HypreParMatrix *
RAP(
const HypreParMatrix *A,
const HypreParMatrix *P);
570 HypreParMatrix *
RAP(
const HypreParMatrix * Rt,
const HypreParMatrix *A,
571 const HypreParMatrix *P);
582 Array2D<double> *blockCoeff=NULL);
587 void EliminateBC(HypreParMatrix &A, HypreParMatrix &Ae,
588 const Array<int> &ess_dof_list,
const Vector &X, Vector &B);
730 virtual operator HYPRE_Solver()
const = 0;
733 virtual HYPRE_PtrToParSolverFcn
SetupFcn()
const = 0;
735 virtual HYPRE_PtrToParSolverFcn
SolveFcn()
const = 0;
738 {
mfem_error(
"HypreSolvers do not support SetOperator!"); }
762 HYPRE_Solver pcg_solver;
795 HYPRE_ParCSRPCGGetNumIterations(pcg_solver, &num_it);
800 virtual operator HYPRE_Solver()
const {
return pcg_solver; }
804 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSetup; }
807 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRPCGSolve; }
820 HYPRE_Solver gmres_solver;
825 void SetDefaultOptions();
850 virtual operator HYPRE_Solver()
const {
return gmres_solver; }
854 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSetup; }
857 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRGMRESSolve; }
870 HYPRE_Solver fgmres_solver;
875 void SetDefaultOptions();
900 virtual operator HYPRE_Solver()
const {
return fgmres_solver; }
904 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRFlexGMRESSetup; }
907 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRFlexGMRESSolve; }
920 virtual operator HYPRE_Solver()
const {
return NULL; }
923 {
return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentitySetup; }
925 {
return (HYPRE_PtrToParSolverFcn) hypre_ParKrylovIdentity; }
936 virtual operator HYPRE_Solver()
const {
return NULL; }
941 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScaleSetup; }
943 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParCSRDiagScale; }
953 HYPRE_Solver sai_precond;
956 void SetDefaultOptions();
961 void ResetSAIPrecond(MPI_Comm comm);
973 virtual operator HYPRE_Solver()
const {
return sai_precond; }
976 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSetup; }
978 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ParaSailsSolve; }
994 HYPRE_Solver euc_precond;
997 void SetDefaultOptions();
1002 void ResetEuclidPrecond(MPI_Comm comm);
1012 virtual operator HYPRE_Solver()
const {
return euc_precond; }
1015 {
return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSetup; }
1017 {
return (HYPRE_PtrToParSolverFcn) HYPRE_EuclidSolve; }
1022 #if MFEM_HYPRE_VERSION >= 21900
1040 HYPRE_Solver ilu_precond;
1043 void SetDefaultOptions();
1049 void ResetILUPrecond();
1064 virtual operator HYPRE_Solver()
const {
return ilu_precond; }
1070 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ILUSetup; }
1074 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ILUSolve; }
1082 HYPRE_Solver amg_precond;
1091 void RecomputeRBMs();
1094 void SetDefaultOptions();
1099 void ResetAMGPrecond();
1120 { HYPRE_BoomerAMGSetPrintLevel(amg_precond, print_level); }
1123 virtual operator HYPRE_Solver()
const {
return amg_precond; }
1126 {
return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSetup; }
1128 {
return (HYPRE_PtrToParSolverFcn) HYPRE_BoomerAMGSolve; }
1134 HypreParMatrix*
DiscreteGrad(ParFiniteElementSpace *edge_fespace,
1135 ParFiniteElementSpace *vert_fespace);
1137 HypreParMatrix*
DiscreteCurl(ParFiniteElementSpace *face_fespace,
1138 ParFiniteElementSpace *edge_fespace);
1169 virtual operator HYPRE_Solver()
const {
return ams; }
1172 {
return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSetup; }
1174 {
return (HYPRE_PtrToParSolverFcn) HYPRE_AMSSolve; }
1209 virtual operator HYPRE_Solver()
const {
return ads; }
1212 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSetup; }
1214 {
return (HYPRE_PtrToParSolverFcn) HYPRE_ADSSolve; }
1254 HYPRE_Solver lobpcg_solver;
1257 mv_InterfaceInterpreter interpreter;
1260 HYPRE_MatvecFunctions matvec_fn;
1266 class HypreMultiVector;
1269 HypreMultiVector * multi_vec;
1278 class HypreMultiVector
1282 mv_MultiVectorPtr mv_ptr;
1292 mv_InterfaceInterpreter & interpreter);
1293 ~HypreMultiVector();
1296 void Randomize(HYPRE_Int seed);
1304 operator mv_MultiVectorPtr()
const {
return mv_ptr; }
1306 mv_MultiVectorPtr & GetMultiVector() {
return mv_ptr; }
1309 static void * OperatorMatvecCreate(
void *A,
void *x );
1310 static HYPRE_Int OperatorMatvec(
void *matvec_data,
1311 HYPRE_Complex
alpha,
1316 static HYPRE_Int OperatorMatvecDestroy(
void *matvec_data );
1318 static HYPRE_Int PrecondSolve(
void *solver,
1322 static HYPRE_Int PrecondSetup(
void *solver,
1390 HYPRE_Solver ame_solver;
1396 HYPRE_Real * eigenvalues;
1399 HYPRE_ParVector * multi_vec;
1403 void createDummyVectors();
1435 #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)
HypreADS(ParFiniteElementSpace *face_fespace)
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.
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
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)
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?
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.
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
Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result a...
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.
Wrapper for Hypre's native parallel ILU preconditioner.
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)
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.
void GetEigenvalues(Array< double > &eigenvalues)
Collect the converged eigenvalues.
Issue warnings on hypre errors.
void SetPreconditioner(HypreSolver &precond)
HypreILU()
Constructor; sets the default options.
HypreAMS(ParFiniteElementSpace *edge_fespace)
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.
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.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
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 SetZeroInitialIterate()
non-hypre setting
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 BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
HypreFGMRES(MPI_Comm comm)
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's FGMRES.
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)
HYPRE_Int NNZ() const
Returns the global number of nonzeros.
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.
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 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 matrix A...
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.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
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
Return the global number of columns.
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)
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 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)
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
PCG Solve function.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
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.
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.
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)
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)
void SetMaxIter(int max_iter)
HypreParMatrix * GetData()
HYPRE_Int * RowPart()
Returns the row partitioning.
void SetPolyOptions(int poly_order, double poly_fraction, int eig_est_cg_iter=10)
Set parameters for polynomial smoothing.
void SetSOROptions(double relax_weight, double omega)
Set SOR-related parameters.
void BooleanMult(int alpha, const int *x, int beta, int *y)
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)
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 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.
A class to initialize the size of a Tensor.
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)
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 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.
void SetPrintLevel(HYPRE_Int print_level)
Set the print level: 0 = none, 1 = setup, 2 = solve, 3 = setup+solve.
HypreParVector(MPI_Comm comm, HYPRE_Int glob_size, HYPRE_Int *col)
Creates vector with given global size and parallel partitioning of the rows/columns given by col...
void SetSubSpaceProjector(Operator &proj)
virtual MFEM_DEPRECATED N_Vector ToNVector()
(DEPRECATED) Return a new wrapper SUNDIALS N_Vector of type SUNDIALS_NVEC_SERIAL. ...
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void SetData(double *_data)
Sets the data of the Vector and the hypre_ParVector to _data.
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.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void SetMassMatrix(HypreParMatrix &M)
HYPRE_Int M() const
Returns the global number of rows.
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.
void GetNumIterations(int &num_iterations)
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function
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...
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.
virtual HYPRE_PtrToParSolverFcn SetupFcn() const
ILU Setup function.
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.
hypre_ParCSRMatrix * StealData()
Changes the ownership of the the matrix.
void SetPrecondUsageMode(int pcg_mode)
void SetSystemsOptions(int dim, bool order_bynodes=false)
virtual ~HypreParMatrix()
Calls hypre's destroy function.
HypreParMatrix * A
The linear system matrix.
ID for class HypreParMatrix.
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
Return the parallel column partitioning array.
HypreParVector & GetEigenvector(unsigned int i)
Extract a single eigenvector.
double omega
SOR parameter (usually in (0,2))
HYPRE_Int * GetRowStarts() const
Return the parallel row partitioning array.
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 SetupFcn() const
hypre's internal Setup function
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 matrix A.
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.
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.
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 parallel row/column partitioning.
double lambda
Taubin's lambda-mu method parameters.
MFEM_DEPRECATED void SetZeroInintialIterate()
deprecated: use SetZeroInitialIterate()
virtual HYPRE_PtrToParSolverFcn SolveFcn() const
hypre's internal Solve function