15 #include "../config/config.hpp"
22 #ifdef MFEM_USE_SUITESPARSE
87 bool final=
false)
const;
119 {
return dot_prod_type == 0 ? MPI_COMM_NULL : comm; }
137 const double damping=1.0);
145 const double damping=1.0);
156 const double damping;
178 int order,
double max_eig_estimate);
190 int order, MPI_Comm comm = MPI_COMM_NULL,
int power_iterations = 10,
191 double power_tolerance = 1e-8);
195 int order,
int power_iterations = 10,
double power_tolerance = 1e-8);
213 double max_eig_estimate;
220 mutable Vector helperVector;
247 void SLI(
const Operator &A,
const Vector &
b, Vector &x,
248 int print_iter = 0,
int max_num_iter = 1000,
249 double RTOLERANCE = 1e-12,
double ATOLERANCE = 1e-24);
252 void SLI(
const Operator &A, Solver &B,
const Vector &
b, Vector &x,
253 int print_iter = 0,
int max_num_iter = 1000,
254 double RTOLERANCE = 1e-12,
double ATOLERANCE = 1e-24);
279 void CG(
const Operator &A,
const Vector &
b, Vector &x,
280 int print_iter = 0,
int max_num_iter = 1000,
281 double RTOLERANCE = 1e-12,
double ATOLERANCE = 1e-24);
284 void PCG(
const Operator &A, Solver &B,
const Vector &
b, Vector &x,
285 int print_iter = 0,
int max_num_iter = 1000,
286 double RTOLERANCE = 1e-12,
double ATOLERANCE = 1e-24);
327 int GMRES(
const Operator &A, Vector &x,
const Vector &
b, Solver &M,
328 int &max_iter,
int m,
double &tol,
double atol,
int printit);
331 void GMRES(
const Operator &A, Solver &B,
const Vector &
b, Vector &x,
332 int print_iter = 0,
int max_num_iter = 1000,
int m = 50,
333 double rtol = 1e-12,
double atol = 1e-24);
358 int BiCGSTAB(
const Operator &A, Vector &x,
const Vector &
b, Solver &M,
359 int &max_iter,
double &tol,
double atol,
int printit);
362 void BiCGSTAB(
const Operator &A, Solver &B,
const Vector &
b, Vector &x,
363 int print_iter = 0,
int max_num_iter = 1000,
364 double rtol = 1e-12,
double atol = 1e-24);
393 void MINRES(
const Operator &A,
const Vector &
b, Vector &x,
int print_it = 0,
394 int max_it = 1000,
double rtol = 1e-12,
double atol = 1e-24);
397 void MINRES(
const Operator &A, Solver &B,
const Vector &
b, Vector &x,
398 int print_it = 0,
int max_it = 1000,
399 double rtol = 1e-12,
double atol = 1e-24);
460 { MFEM_WARNING(
"L-BFGS won't use the given preconditioner."); }
462 { MFEM_WARNING(
"L-BFGS won't use the given solver."); }
469 int aGMRES(
const Operator &A, Vector &x,
const Vector &
b,
470 const Operator &M,
int &max_iter,
471 int m_max,
int m_min,
int m_step,
double cf,
472 double &tol,
double &atol,
int printit);
508 { MFEM_ABORT(
"The objective gradient is not implemented."); }
546 { MFEM_ABORT(
"Not meaningful for this solver."); }
548 { MFEM_ABORT(
"Not meaningful for this solver."); }
724 bool final)
override;
728 #ifdef MFEM_USE_SUITESPARSE
743 mutable double Info[UMFPACK_INFO];
800 #endif // MFEM_USE_SUITESPARSE
804 #endif // MFEM_SOLVERS
const Array< int > * ess_dofs_list
Not owned.
const Vector * GetInequalityVec_Hi() const
Conjugate gradient method.
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
double Info[UMFPACK_INFO]
void Mult(const Vector &b, Vector &x) const
Solve the system LUx = b, where L and U are the block ILU factors.
Chebyshev accelerated smoothing with given vector, no matrix necessary.
int GetNumIterations() const
MPI_Comm GetComm() const
Return the associated MPI communicator, or MPI_COMM_NULL if no communicator is set.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
void SetHistorySize(int dim)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void SetSize(int s)
Resize the vector to size s.
void Monitor(int it, double norm, const Vector &r, const Vector &x, bool final=false) const
virtual void ProcessNewState(const Vector &x) const
This method can be overloaded in derived classes to perform computations that need knowledge of the n...
NewtonSolver(MPI_Comm _comm)
~OperatorChebyshevSmoother()
class IterativeSolver * iter_solver
The last IterativeSolver to which this monitor was attached.
T * GetData()
Returns the data.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Direct sparse solver using KLU.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
virtual void MultTranspose(const Vector &b, Vector &x) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
UMFPackSolver(SparseMatrix &A, bool _use_long_ints=false)
Factorize the given SparseMatrix using the defaults. For larger matrices, if the solver fails...
void SetOperator(const Operator &op)
const Operator * GetC() const
Reordering
The reordering method used by the BlockILU factorization.
virtual ~IterativeSolverMonitor()
ResidualBCMonitor(const Array< int > &ess_dofs_list_)
virtual double CalcObjective(const Vector &x) const =0
Objective F(x). In parallel, the result should be reduced over tasks.
void add(const Vector &v1, const Vector &v2, Vector &v)
virtual ~OptimizationSolver()
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 ...
double Norm(const Vector &x) const
double GetFinalNorm() const
int BiCGSTAB(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, double &tol, double atol, int printit)
BiCGSTAB method. (tolerances are squared)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
FGMRESSolver(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 Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
void SetLinearConstraint(const Vector &_w, double _a)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
void SetPrintLevel(int print_lvl)
const Vector * GetBoundsVec_Lo() const
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
double solve(double l, const Vector &xt, Vector &x, int &nclip) const
Solve QP at fixed lambda.
int GetNumConstraints() const
OperatorChebyshevSmoother(Operator *oper_, const Vector &d, const Array< int > &ess_tdof_list, int order, double max_eig_estimate)
void MINRES(const Operator &A, const Vector &b, Vector &x, int print_it, int max_it, double rtol, double atol)
MINRES method without preconditioner. (tolerances are squared)
LBFGSSolver(MPI_Comm _comm)
Jacobi smoothing for a given bilinear form (no matrix necessary).
~OperatorJacobiSmoother()
void SetMaxIter(int max_it)
void median(const Vector &lo, const Vector &hi)
v = median(v,lo,hi) entrywise. Implementation assumes lo <= hi.
void SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
Newton's method for solving F(x)=b for a given operator F.
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
void print_iteration(int it, double r, double l) const
void Setup(const Vector &diag)
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
double Dot(const Vector &x, const Vector &y) const
BiCGSTABSolver(MPI_Comm _comm)
void SetIterativeSolver(const IterativeSolver &solver)
This method is invoked by ItertiveSolver::SetMonitor, informing the monitor which IterativeSolver is ...
const Vector * GetBoundsVec_Hi() const
SLBQPOptimizer(MPI_Comm _comm)
BlockILU(int block_size_, Reordering reordering_=Reordering::MINIMUM_DISCARDED_FILL, int k_fill_=0)
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
OperatorJacobiSmoother(const BilinearForm &a, const Array< int > &ess_tdof_list, const double damping=1.0)
void SetPrintLevel(int print_lvl)
Set the print level field in the Control data member.
Stationary linear iteration: x <- x + B (b - A x)
void SetInequalityConstraint(const Vector &dl, const Vector &dh)
void SetEqualityConstraint(const Vector &c)
Abstract base class for an iterative solver monitor.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
void SetOperator(const Operator &op_)
Set/update the solver for the given operator.
void SetAbsTol(double atol)
void SetRelTol(double rtol)
double Control[UMFPACK_CONTROL]
Monitor that checks whether the residual is zero at a given set of dofs.
Abstract base class for iterative solver.
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
OptimizationSolver(MPI_Comm _comm)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
SLISolver(MPI_Comm _comm)
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
GMRESSolver(MPI_Comm _comm)
const Vector * GetEqualityVec() const
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
void SetMonitor(IterativeSolverMonitor &m)
Set the iterative solver monitor.
void SetBounds(const Vector &_lo, const Vector &_hi)
const Operator * GetD() const
Abstract solver for OptimizationProblems.
virtual double ComputeScalingFactor(const Vector &x, const Vector &b) const
This method can be overloaded in derived classes to implement line search algorithms.
OptimizationProblem(int insize, const Operator *C_, const Operator *D_)
In parallel, insize is the number of the local true dofs.
void SetSolutionBounds(const Vector &xl, const Vector &xh)
const Vector * GetInequalityVec_Lo() const
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual void Mult(const Vector &xt, Vector &x) const =0
Operator application: y=A(x).
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void Mult(const Vector &xt, Vector &x) const
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 void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void CalcObjectiveGrad(const Vector &x, Vector &grad) const
The result grad is expected to enter with the correct size.
void MonitorResidual(int it, double norm, const Vector &r, bool final) override
Monitor the residual vector r.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
MINRESSolver(MPI_Comm _comm)
virtual void MonitorSolution(int it, double norm, const Vector &x, bool final)
Monitor the solution vector x.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
KLUSolver(SparseMatrix &A)
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, double cf, double &tol, double &atol, int printit)
const OptimizationProblem * problem
Rank 3 tensor (array of matrices)
virtual void SetOptimizationProblem(const OptimizationProblem &prob)
IterativeSolverMonitor * monitor
virtual void MonitorResidual(int it, double norm, const Vector &r, bool final)
Monitor the residual vector r.
int width
Dimension of the input / number of columns in the matrix.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
const Operator * C
Not owned, some can remain unused (NULL).
void SLI(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Stationary linear iteration. (tolerances are squared)
int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, int m, double &tol, double atol, int printit)
GMRES method. (tolerances are squared)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.