|
ROL
|
Provides an interface to run the trust-region algorithm of Kelley and Sachs. More...
#include <ROL_TypeB_KelleySachsAlgorithm.hpp>
Inheritance diagram for ROL::TypeB::KelleySachsAlgorithm< Real >:Public Member Functions | |
| KelleySachsAlgorithm (ParameterList &list, const Ptr< Secant< Real > > &secant=nullPtr) | |
| void | run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout) override |
| Run algorithm on bound constrained problems (Type-B). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method. | |
| void | run (Problem< Real > &problem, std::ostream &outStream=std::cout) override |
| Run algorithm on bound constrained problems (Type-B). This is the primary Type-B interface. | |
| void | run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, const Vector< Real > &linear_eres, std::ostream &outStream=std::cout) override |
| Run algorithm on bound constrained problems with explicit linear constraints (Type-B). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method. | |
| void | writeHeader (std::ostream &os) const override |
| Print iterate header. | |
| void | writeName (std::ostream &os) const override |
| Print step name. | |
| void | writeOutput (std::ostream &os, const bool write_header=false) const override |
| Print iterate status. | |
Public Member Functions inherited from ROL::TypeB::Algorithm< Real > | |
| virtual | ~Algorithm () |
| Algorithm () | |
| Constructor, given a step and a status test. | |
| void | setStatusTest (const Ptr< StatusTest< Real > > &status, const bool combineStatus=false) |
| virtual void | run (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout) |
| Run algorithm on bound constrained problems (Type-B). This is the primary Type-B interface. | |
| virtual void | run (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, std::ostream &outStream=std::cout) |
| Run algorithm on bound constrained problems with explicit linear constraints (Type-B). This is the primary Type-B with explicit linear constraints interface. | |
| virtual void | run (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &linear_icon, Vector< Real > &linear_imul, BoundConstraint< Real > &linear_ibnd, std::ostream &outStream=std::cout) |
| Run algorithm on bound constrained problems with explicit linear constraints (Type-B). This is the primary Type-B with explicit linear constraints interface. | |
| virtual void | run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &linear_icon, Vector< Real > &linear_imul, BoundConstraint< Real > &linear_ibnd, const Vector< Real > &linear_ires, std::ostream &outStream=std::cout) |
| Run algorithm on bound constrained problems with explicit linear constraints (Type-B). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method. | |
| virtual void | run (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, Constraint< Real > &linear_icon, Vector< Real > &linear_imul, BoundConstraint< Real > &linear_ibnd, std::ostream &outStream=std::cout) |
| Run algorithm on bound constrained problems with explicit linear constraints (Type-B). This is the primary Type-B with explicit linear constraints interface. | |
| virtual void | run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, const Vector< Real > &linear_eres, Constraint< Real > &linear_icon, Vector< Real > &linear_imul, BoundConstraint< Real > &linear_ibnd, const Vector< Real > &linear_ires, std::ostream &outStream=std::cout) |
| Run algorithm on bound constrained problems with explicit linear constraints (Type-B). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method. | |
| virtual void | writeExitStatus (std::ostream &os) const |
| Ptr< const AlgorithmState< Real > > | getState () const |
| void | reset () |
Private Member Functions | |
| void | initialize (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout) |
| Real | trqsol (const Real xtx, const Real ptp, const Real ptx, const Real del) const |
| Real | trpcg (Vector< Real > &w, int &iflag, int &iter, const Vector< Real > &g, const Vector< Real > &x, const Vector< Real > &g0, const Real del, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, Real eps, const Real tol, const int itermax, Vector< Real > &p, Vector< Real > &q, Vector< Real > &r, Vector< Real > &t, Vector< Real > &pwa, Vector< Real > &dwa, std::ostream &outStream=std::cout) const |
| void | applyFreeHessian (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real eps, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, Real &tol, Vector< Real > &pwa, Vector< Real > &dwa) const |
| void | applyFreePrecond (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real eps, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, Real &tol, Vector< Real > &pwa, Vector< Real > &dwa) const |
Private Attributes | |
| Ptr< TrustRegionModel_U< Real > > | model_ |
| Container for trust-region model. | |
| Real | delMax_ |
| Maximum trust-region radius (default: ROL_INF) | |
| Real | eta0_ |
| Step acceptance threshold (default: 0.05) | |
| Real | eta1_ |
| Radius decrease threshold (default: 0.05) | |
| Real | eta2_ |
| Radius increase threshold (default: 0.9) | |
| Real | gamma0_ |
| Radius decrease rate (negative rho) (default: 0.0625) | |
| Real | gamma1_ |
| Radius decrease rate (positive rho) (default: 0.25) | |
| Real | gamma2_ |
| Radius increase rate (default: 2.5) | |
| Real | TRsafe_ |
| Safeguard size for numerically evaluating ratio (default: 1e2) | |
| Real | eps_ |
| Safeguard for numerically evaluating ratio. | |
| TRUtils::ETRFlag | TRflag_ |
| Trust-region exit flag. | |
| int | SPflag_ |
| Subproblem solver termination flag. | |
| int | SPiter_ |
| Subproblem solver iteration count. | |
| ESecant | esec_ |
| Secant type (default: Limited-Memory BFGS) | |
| bool | useSecantPrecond_ |
| Flag to use secant as a preconditioner (default: false) | |
| bool | useSecantHessVec_ |
| Flag to use secant as Hessian (default: false) | |
| int | maxit_ |
| Maximum number of CG iterations (default: 20) | |
| Real | tol1_ |
| Absolute tolerance for truncated CG (default: 1e-4) | |
| Real | tol2_ |
| Relative tolerance for truncated CG (default: 1e-2) | |
| int | minit_ |
| Maximum number of minor (subproblem solve) iterations (default: 10) | |
| Real | mu0_ |
| Sufficient decrease parameter (default: 1e-2) | |
| Real | mu1_ |
| Sufficient decrease parameter postsmoothing (default: 0.9999) | |
| Real | eps0_ |
| Epsilon binding set tolerance (default: 1e-3) | |
| Real | beta_ |
| Projected search backtracking rate (default: 1e-2) | |
| Real | alpha0_ |
| Initial step size for projected search (default: 1) | |
| int | nhess_ |
| Number of Hessian applications. | |
| unsigned | verbosity_ |
| Output level (default: 0) | |
| bool | writeHeader_ |
| Flag to write header at every iteration. | |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::TypeB::Algorithm< Real > | |
| void | initialize (const Vector< Real > &x, const Vector< Real > &g) |
| Real | optimalityCriterion (const Vector< Real > &x, const Vector< Real > &g, Vector< Real > &primal, std::ostream &outStream=std::cout) const |
Protected Attributes inherited from ROL::TypeB::Algorithm< Real > | |
| const Ptr< CombinedStatusTest< Real > > | status_ |
| const Ptr< AlgorithmState< Real > > | state_ |
| Ptr< PolyhedralProjection< Real > > | proj_ |
Provides an interface to run the trust-region algorithm of Kelley and Sachs.
Definition at line 27 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
| ROL::TypeB::KelleySachsAlgorithm< Real >::KelleySachsAlgorithm | ( | ParameterList & | list, |
| const Ptr< Secant< Real > > & | secant = nullPtr |
||
| ) |
Definition at line 17 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
References ROL::StringToESecant().
|
overridevirtual |
Run algorithm on bound constrained problems (Type-B). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method.
Implements ROL::TypeB::Algorithm< Real >.
Definition at line 94 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
References ROL::Accept, ROL::Vector< Real >::clone(), ROL::Objective< Real >::gradient(), ROL::TRUtils::NPOSPREDNEG, ROL::Vector< Real >::plus(), ROL::TRUtils::POSPREDNEG, ROL::BoundConstraint< Real >::project(), ROL::BoundConstraint< Real >::pruneActive(), ROL::TRUtils::QMINSUFDEC, ROL::Revert, ROL::Vector< Real >::set(), ROL::TRUtils::SUCCESS, ROL::Trial, ROL::TRUtils::TRNAN, ROL::Objective< Real >::update(), ROL::Objective< Real >::value(), and ROL::TypeB::Algorithm< Real >::writeExitStatus().
|
overridevirtual |
Run algorithm on bound constrained problems (Type-B). This is the primary Type-B interface.
Reimplemented from ROL::TypeB::Algorithm< Real >.
Definition at line 258 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
References ROL::Problem< Real >::getPolyhedralProjection(), and ROL::TypeB::Algorithm< Real >::run().
|
overridevirtual |
Run algorithm on bound constrained problems with explicit linear constraints (Type-B). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method.
Reimplemented from ROL::TypeB::Algorithm< Real >.
Definition at line 269 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
|
overridevirtual |
Print iterate header.
Reimplemented from ROL::TypeB::Algorithm< Real >.
Definition at line 455 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
References ROL::CG_FLAG_SUCCESS, ROL::CG_FLAG_UNDEFINED, ROL::ECGFlagToString(), ROL::TRUtils::ETRFlagToString(), ROL::NumberToString(), ROL::TRUtils::SUCCESS, and ROL::TRUtils::UNDEFINED.
|
overridevirtual |
Print step name.
Reimplemented from ROL::TypeB::Algorithm< Real >.
Definition at line 500 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
|
overridevirtual |
Print iterate status.
Reimplemented from ROL::TypeB::Algorithm< Real >.
Definition at line 507 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
|
private |
|
private |
Definition at line 281 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
References zero.
|
private |
Definition at line 303 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
References ROL::Vector< Real >::apply(), ROL::Vector< Real >::axpy(), ROL::Vector< Real >::norm(), ROL::Vector< Real >::plus(), ROL::Vector< Real >::scale(), ROL::Vector< Real >::set(), zero, and ROL::Vector< Real >::zero().
|
private |
Definition at line 392 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
References ROL::Vector< Real >::dual(), ROL::TrustRegionModel_U< Real >::hessVec(), ROL::Vector< Real >::plus(), ROL::BoundConstraint< Real >::pruneActive(), ROL::BoundConstraint< Real >::pruneInactive(), and ROL::Vector< Real >::set().
|
private |
Definition at line 424 of file ROL_TypeB_KelleySachsAlgorithm_Def.hpp.
References ROL::Vector< Real >::dual(), ROL::Vector< Real >::plus(), ROL::TrustRegionModel_U< Real >::precond(), ROL::BoundConstraint< Real >::pruneActive(), ROL::BoundConstraint< Real >::pruneInactive(), and ROL::Vector< Real >::set().
|
private |
Container for trust-region model.
Definition at line 29 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Maximum trust-region radius (default: ROL_INF)
Definition at line 32 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Step acceptance threshold (default: 0.05)
Definition at line 33 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Radius decrease threshold (default: 0.05)
Definition at line 34 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Radius increase threshold (default: 0.9)
Definition at line 35 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Radius decrease rate (negative rho) (default: 0.0625)
Definition at line 36 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Radius decrease rate (positive rho) (default: 0.25)
Definition at line 37 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Radius increase rate (default: 2.5)
Definition at line 38 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Safeguard size for numerically evaluating ratio (default: 1e2)
Definition at line 39 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Safeguard for numerically evaluating ratio.
Definition at line 40 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Trust-region exit flag.
Definition at line 43 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Subproblem solver termination flag.
Definition at line 44 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Subproblem solver iteration count.
Definition at line 45 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Secant type (default: Limited-Memory BFGS)
Definition at line 48 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Flag to use secant as a preconditioner (default: false)
Definition at line 49 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Flag to use secant as Hessian (default: false)
Definition at line 50 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Maximum number of CG iterations (default: 20)
Definition at line 53 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Absolute tolerance for truncated CG (default: 1e-4)
Definition at line 54 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Relative tolerance for truncated CG (default: 1e-2)
Definition at line 55 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Maximum number of minor (subproblem solve) iterations (default: 10)
Definition at line 58 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Sufficient decrease parameter (default: 1e-2)
Definition at line 59 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Sufficient decrease parameter postsmoothing (default: 0.9999)
Definition at line 60 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Epsilon binding set tolerance (default: 1e-3)
Definition at line 61 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Projected search backtracking rate (default: 1e-2)
Definition at line 62 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Initial step size for projected search (default: 1)
Definition at line 63 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
mutableprivate |
Number of Hessian applications.
Definition at line 65 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Output level (default: 0)
Definition at line 66 of file ROL_TypeB_KelleySachsAlgorithm.hpp.
|
private |
Flag to write header at every iteration.
Definition at line 67 of file ROL_TypeB_KelleySachsAlgorithm.hpp.