24#ifdef MFEM_USE_SUITESPARSE
117 MPI_Comm comm = MPI_COMM_NULL;
184 bool final=
false)
const;
303 {
return dot_prod_type == 0 ? MPI_COMM_NULL : comm; }
331 const real_t damping=1.0);
344 const real_t damping=1.0);
377 bool use_abs_diag =
false;
383 const bool allow_updates;
404 int order,
real_t max_eig_estimate);
410 int order,
real_t max_eig_estimate);
422 int order, MPI_Comm comm = MPI_COMM_NULL,
423 int power_iterations = 10,
424 real_t power_tolerance = 1e-8);
430 int order, MPI_Comm comm = MPI_COMM_NULL,
431 int power_iterations = 10,
432 real_t power_tolerance = 1e-8);
436 int order,
int power_iterations = 10,
437 real_t power_tolerance = 1e-8);
443 int order,
int power_iterations = 10,
444 real_t power_tolerance = 1e-8);
473 mutable Vector helperVector;
501void SLI(
const Operator &A,
const Vector &
b, Vector &x,
502 int print_iter = 0,
int max_num_iter = 1000,
506void SLI(
const Operator &A, Solver &B,
const Vector &
b, Vector &x,
507 int print_iter = 0,
int max_num_iter = 1000,
535void CG(
const Operator &A,
const Vector &
b, Vector &x,
536 int print_iter = 0,
int max_num_iter = 1000,
540void PCG(
const Operator &A, Solver &B,
const Vector &
b, Vector &x,
541 int print_iter = 0,
int max_num_iter = 1000,
585int GMRES(
const Operator &A, Vector &x,
const Vector &
b, Solver &M,
586 int &max_iter,
int m,
real_t &tol,
real_t atol,
int printit);
589void GMRES(
const Operator &A, Solver &B,
const Vector &
b, Vector &x,
590 int print_iter = 0,
int max_num_iter = 1000,
int m = 50,
617int BiCGSTAB(
const Operator &A, Vector &x,
const Vector &
b, Solver &M,
621void BiCGSTAB(
const Operator &A, Solver &B,
const Vector &
b, Vector &x,
622 int print_iter = 0,
int max_num_iter = 1000,
653void MINRES(
const Operator &A,
const Vector &
b, Vector &x,
int print_it = 0,
654 int max_it = 1000,
real_t rtol = 1e-12,
real_t atol = 1e-24);
657void MINRES(
const Operator &A, Solver &B,
const Vector &
b, Vector &x,
658 int print_it = 0,
int max_it = 1000,
695 const real_t fnorm)
const;
702 const real_t fnorm)
const;
744 const real_t rtol_max = 0.9,
759 for (
int i = 0; i <
skArray.Size(); i++)
771 for (
int i = 0; i <
m; i++)
804 { MFEM_WARNING(
"L-BFGS won't use the given preconditioner."); }
806 { MFEM_WARNING(
"L-BFGS won't use the given solver."); }
816int aGMRES(
const Operator &A, Vector &x,
const Vector &
b,
817 const Operator &M,
int &max_iter,
818 int m_max,
int m_min,
int m_step,
real_t cf,
822class HiopOptimizationProblem;
849 mutable bool new_x =
true;
860 bool NewX()
const {
return new_x; }
872 { MFEM_ABORT(
"The objective gradient is not implemented."); }
910 { MFEM_ABORT(
"Not meaningful for this solver."); }
912 { MFEM_ABORT(
"Not meaningful for this solver."); }
1088 bool final)
override;
1092#ifdef MFEM_USE_SUITESPARSE
1179 std::unique_ptr<DenseMatrixInverse[]> block_solvers;
1199 bool ownA,
bool ownS0,
bool ownS1)
1200 :
Solver(A_->
NumRows()), A(A_, ownA), S0(S0_, ownS0), S1(S1_, ownS1) { }
1223 const bool parallel;
1225 mutable int global_size;
1254 Solver *solver =
nullptr;
1258 void Orthogonalize(
const Vector &v,
Vector &v_ortho)
const;
1273 void Mult(
const Vector &x,
Vector &y,
bool transpose)
const;
1276 bool op_is_symmetric =
true,
bool own_aux_map =
false);
1279 {
Mult(x, y,
true); }
1286#ifdef MFEM_USE_LAPACK
1324 { res_change_termination_tol_ = tol; }
1388 mutable int max_nnz_;
1395 real_t res_change_termination_tol_;
1405 mutable bool NNLS_qrres_on_;
1408 mutable Vector row_scaling_;
T * GetData()
Returns the data.
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 ...
HypreSmoother & GetSmoother()
AuxSpaceSmoother(const HypreParMatrix &op, HypreParMatrix *aux_map, bool op_is_symmetric=true, bool own_aux_map=false)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the BiCGSTAB method.
BiCGSTABSolver(MPI_Comm comm_)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Reordering
The reordering method used by the BlockILU factorization.
BlockILU(int block_size_, Reordering reordering_=Reordering::MINIMUM_DISCARDED_FILL, int k_fill_=0)
void Mult(const Vector &b, Vector &x) const
Solve the system LUx = b, where L and U are the block ILU factors.
void SetOperator(const Operator &op)
Conjugate gradient method.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Data type dense matrix using column-major storage.
Rank 3 tensor (array of matrices)
Block diagonal solver for A, each block is inverted by direct solver.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
DirectSubBlockSolver(const SparseMatrix &A, const SparseMatrix &block_dof)
virtual void Mult(const Vector &x, Vector &y) const
Direct solution of the block diagonal linear system.
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the FGMRES method.
FGMRESSolver(MPI_Comm comm_)
void SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the GMRES method.
GMRESSolver(MPI_Comm comm_)
Internal class - adapts the OptimizationProblem class to HiOp's interface.
Wrapper for hypre's ParCSR matrix class.
Parallel smoothers in hypre.
Abstract base class for an iterative solver monitor.
void SetIterativeSolver(const IterativeSolver &solver)
This method is invoked by IterativeSolver::SetMonitor, informing the monitor which IterativeSolver is...
const class IterativeSolver * iter_solver
The last IterativeSolver to which this monitor was attached.
virtual void MonitorSolution(int it, real_t norm, const Vector &x, bool final)
Monitor the solution vector x.
virtual ~IterativeSolverMonitor()
virtual void MonitorResidual(int it, real_t norm, const Vector &r, bool final)
Monitor the residual vector r.
Abstract base class for iterative solver.
real_t abs_tol
Absolute tolerance.
void Monitor(int it, real_t norm, const Vector &r, const Vector &x, bool final=false) const
Monitor both the residual r and the solution x.
real_t rel_tol
Relative tolerance.
PrintLevel print_options
Output behavior for the iterative solver.
real_t GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult().
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
virtual real_t Dot(const Vector &x, const Vector &y) const
Return the standard (l2, i.e., Euclidean) inner product of x and y.
void SetRelTol(real_t rtol)
real_t GetFinalRelNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult(),...
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
int print_level
(DEPRECATED) Legacy print level definition, which is left for compatibility with custom iterative sol...
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
int max_iter
Limit for the number of iterations the solver is allowed to do.
void SetMaxIter(int max_it)
MPI_Comm GetComm() const
Return the associated MPI communicator, or MPI_COMM_NULL if no communicator is set.
real_t GetInitialNorm() const
Returns the initial residual norm from the last call to Mult().
static int GuessLegacyPrintLevel(PrintLevel)
Use some heuristics to guess a legacy print level corresponding to the given PrintLevel.
IterativeSolverMonitor * monitor
void SetMonitor(IterativeSolverMonitor &m)
Set the iterative solver monitor.
bool GetConverged() const
Returns true if the last call to Mult() converged successfully.
PrintLevel FromLegacyPrintLevel(int)
Convert a legacy print level integer to a PrintLevel object.
void SetAbsTol(real_t atol)
real_t Norm(const Vector &x) const
Return the inner product norm of x, using the inner product defined by Dot()
Direct sparse solver using KLU.
virtual void MultTranspose(const Vector &b, Vector &x) const
Direct solution of the transposed linear system using KLU.
virtual void Mult(const Vector &b, Vector &x) const
Direct solution of the linear system using KLU.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
KLUSolver(SparseMatrix &A)
LBFGSSolver(MPI_Comm comm_)
Array< Vector * > skArray
void InitializeStorageVectors()
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
void SetHistorySize(int dim)
void DeleteStorageVectors()
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Array< Vector * > ykArray
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
MINRESSolver(MPI_Comm comm_)
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the MINRES method.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
void SetNormalize(bool n)
Set a flag to determine whether to call NormalizeConstraints().
void SetResidualChangeTolerance(real_t tol)
Set threshold on relative change in residual over nStallCheck_ iterations.
void SetRHSDelta(real_t d)
Set RHS vector constant shift, defining rhs_lb and rhs_ub in Solve().
void Solve(const Vector &rhs_lb, const Vector &rhs_ub, Vector &soln) const
Solve the NNLS problem. Specifically, we find a vector soln, such that rhs_lb < mat*soln < rhs_ub is ...
QRresidualMode
Enumerated types of QRresidual mode. Options are 'off': the residual is calculated normally,...
void SetZeroTolerance(real_t tol)
Set the magnitude of projected residual entries that are considered zero. Increasing this value relax...
void SetInnerIterations(int n)
Set the maximum number of inner iterations in Solve().
void Mult(const Vector &w, Vector &sol) const override
Compute the non-negative least squares solution to the underdetermined system.
void SetQRResidualMode(const QRresidualMode qr_residual_mode)
Set the residual calculation mode for the NNLS solver. See QRresidualMode enum above for details.
void SetOperator(const Operator &op) override
The operator must be a DenseMatrix.
void SetStallCheck(int n)
Set the number of iterations to use for stall checking.
void NormalizeConstraints(Vector &rhs_lb, Vector &rhs_ub) const
Normalize the constraints such that the tolerances for each constraint (i.e. (UB - LB)/2) are equal....
void SetTolerance(real_t tol)
Set the target absolute residual norm tolerance for convergence.
void SetOuterIterations(int n)
Set the maximum number of outer iterations in Solve().
void SetMinNNZ(int min_nnz)
Set the minimum number of nonzeros required for the solution.
void SetMaxNNZ(int max_nnz)
Set the maximum number of nonzeros required for the solution, as an early termination condition.
void SetVerbosity(int v)
Set verbosity. If set to 0: print nothing; if 1: just print results; if 2: print short update on ever...
Newton's method for solving F(x)=b for a given operator F.
void SetAdaptiveLinRtol(const int type=2, const real_t rtol0=0.5, const real_t rtol_max=0.9, const real_t alpha=0.5 *(1.0+sqrt(5.0)), const real_t gamma=1.0)
Enable adaptive linear solver relative tolerance algorithm.
void AdaptiveLinRtolPreSolve(const Vector &x, const int it, const real_t fnorm) const
Method for the adaptive linear solver rtol invoked before the linear solve.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
NewtonSolver(MPI_Comm comm_)
virtual real_t ComputeScalingFactor(const Vector &x, const Vector &b) const
This method can be overloaded in derived classes to implement line search algorithms.
void AdaptiveLinRtolPostSolve(const Vector &x, const Vector &b, const int it, const real_t fnorm) const
Method for the adaptive linear solver rtol invoked after the linear solve.
virtual void ProcessNewState(const Vector &x) const
This method can be overloaded in derived classes to perform computations that need knowledge of the n...
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Chebyshev accelerated smoothing with given vector, no matrix necessary.
void SetOperator(const Operator &op_)
Set/update the solver for the given operator.
void Mult(const Vector &x, Vector &y) const
Approach the solution of the linear system by applying Chebyshev smoothing.
OperatorChebyshevSmoother(const Operator &oper_, const Vector &d, const Array< int > &ess_tdof_list, int order, int power_iterations=10, real_t power_tolerance=1e-8)
void MultTranspose(const Vector &x, Vector &y) const
Approach the solution of the transposed linear system by applying Chebyshev smoothing.
OperatorChebyshevSmoother(const Operator &oper_, const Vector &d, const Array< int > &ess_tdof_list, int order, real_t max_eig_estimate)
~OperatorChebyshevSmoother()
Pointer to an Operator of a specified type.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
Jacobi smoothing for a given bilinear form (no matrix necessary).
void Setup(const Vector &diag)
~OperatorJacobiSmoother()
OperatorJacobiSmoother(const real_t damping=1.0)
Default constructor: the diagonal will be computed by subsequent calls to SetOperator() using the Ope...
void SetPositiveDiagonal(bool pos_diag=true)
Replace diagonal entries with their absolute values.
void SetOperator(const Operator &op)
Recompute the diagonal using the method AssembleDiagonal of the given new Operator,...
void MultTranspose(const Vector &x, Vector &y) const
Approach the solution of the transposed linear system by applying Jacobi smoothing.
void Mult(const Vector &x, Vector &y) const
Approach the solution of the linear system by applying Jacobi smoothing.
int width
Dimension of the input / number of columns in the matrix.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
int NumRows() const
Get the number of rows (size of output) of the Operator. Synonym with Height().
virtual void CalcObjectiveGrad(const Vector &x, Vector &grad) const
The result grad is expected to enter with the correct size.
const Vector * GetEqualityVec() const
const Operator * C
Not owned, some can remain unused (NULL).
void SetSolutionBounds(const Vector &xl, const Vector &xh)
void SetEqualityConstraint(const Vector &c)
int GetNumConstraints() const
virtual real_t CalcObjective(const Vector &x) const =0
Objective F(x). In parallel, the result should be reduced over tasks.
const Vector * GetInequalityVec_Hi() const
const Operator * GetC() const
const Vector * GetBoundsVec_Hi() const
const Operator * GetD() const
OptimizationProblem(int insize, const Operator *C_, const Operator *D_)
In parallel, insize is the number of the local true dofs.
const Vector * GetInequalityVec_Lo() const
const Vector * GetBoundsVec_Lo() const
void SetInequalityConstraint(const Vector &dl, const Vector &dh)
Abstract solver for OptimizationProblems.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
virtual ~OptimizationSolver()
virtual void Mult(const Vector &xt, Vector &x) const =0
Operator application: y=A(x).
const OptimizationProblem * problem
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
OptimizationSolver(MPI_Comm comm_)
Solver wrapper which orthogonalizes the input and output vector.
void Mult(const Vector &b, Vector &x) const
Perform the action of the OrthoSolver: P * solver * P where P is the projection to the subspace of ve...
OrthoSolver(MPI_Comm mycomm_)
void SetSolver(Solver &s)
Set the solver used by the OrthoSolver.
virtual void SetOperator(const Operator &op)
Set the Operator that is the OrthoSolver is to invert (approximately).
ProductSolver(Operator *A_, Solver *S0_, Solver *S1_, bool ownA, bool ownS0, bool ownS1)
virtual void MultTranspose(const Vector &x, Vector &y) const
Solution of the transposed linear system using a product of subsolvers.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
virtual void Mult(const Vector &x, Vector &y) const
Solution of the linear system using a product of subsolvers.
Monitor that checks whether the residual is zero at a given set of dofs.
void MonitorResidual(int it, real_t norm, const Vector &r, bool final) override
Monitor the residual vector r.
const Array< int > * ess_dofs_list
Not owned.
ResidualBCMonitor(const Array< int > &ess_dofs_list_)
virtual void Mult(const Vector &xt, Vector &x) const
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
SLBQPOptimizer(MPI_Comm comm_)
void SetBounds(const Vector &lo_, const Vector &hi_)
void SetLinearConstraint(const Vector &w_, real_t a_)
void print_iteration(int it, real_t r, real_t l) const
real_t solve(real_t l, const Vector &xt, Vector &x, int &nclip) const
Solve QP at fixed lambda.
Stationary linear iteration: x <- x + B (b - A x)
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using Stationary Linear Iteration.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
SLISolver(MPI_Comm comm_)
Direct sparse solver using UMFPACK.
UMFPackSolver(bool use_long_ints_=false)
For larger matrices, if the solver fails, set the parameter use_long_ints_ = true.
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
real_t Control[UMFPACK_CONTROL]
virtual void Mult(const Vector &b, Vector &x) const
Direct solution of the linear system using UMFPACK.
virtual void MultTranspose(const Vector &b, Vector &x) const
Direct solution of the transposed linear system using UMFPACK.
real_t Info[UMFPACK_INFO]
void SetPrintLevel(int print_lvl)
Set the print level field in the Control data member.
UMFPackSolver(SparseMatrix &A, bool use_long_ints_=false)
Factorize the given SparseMatrix using the defaults. For larger matrices, if the solver fails,...
void median(const Vector &lo, const Vector &hi)
v = median(v,lo,hi) entrywise. Implementation assumes lo <= hi.
void SetSize(int s)
Resize the vector to size s.
int BiCGSTAB(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, real_t &tol, real_t atol, int printit)
BiCGSTAB method. (tolerances are squared)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
void add(const Vector &v1, const Vector &v2, Vector &v)
int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, int m, real_t &tol, real_t atol, int printit)
GMRES method. (tolerances are squared)
int aGMRES(const Operator &A, Vector &x, const Vector &b, const Operator &M, int &max_iter, int m_max, int m_min, int m_step, real_t cf, real_t &tol, real_t &atol, int printit)
void SLI(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Stationary linear iteration. (tolerances are squared)
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
void MINRES(const Operator &A, const Vector &b, Vector &x, int print_it, int max_it, real_t rtol, real_t atol)
MINRES method without preconditioner. (tolerances are squared)
Settings for the output behavior of the IterativeSolver.
PrintLevel & Iterations()
bool errors
If a fatal problem has been detected the failure will be reported to mfem::err.
bool iterations
Detailed information about each iteration will be reported to mfem::out.
PrintLevel & FirstAndLast()
PrintLevel()=default
Initializes the print level to suppress.
bool warnings
If a non-fatal problem has been detected some context-specific information will be reported to mfem::...
bool first_and_last
Information about the first and last iteration will be printed to mfem::out.
bool summary
A summary of the solver process will be reported after the last iteration to mfem::out.