ROL
ROL_TypeB_ColemanLiAlgorithm.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_TYPEB_COLEMANLIALGORITHM_HPP
11#define ROL_TYPEB_COLEMANLIALGORITHM_HPP
12
18
23namespace ROL {
24namespace TypeB {
25
26template<typename Real>
28private:
29 Ptr<TrustRegionModel_U<Real>> model_;
30
31 // TRUST REGION PARAMETERS
32 Real delMax_;
33 Real eta0_;
34 Real eta1_;
35 Real eta2_;
36 Real gamma0_;
37 Real gamma1_;
38 Real gamma2_;
39 Real TRsafe_;
40 Real eps_;
42
43 // NONMONOTONE PARAMETER
44 bool useNM_;
46
47 // ITERATION FLAGS/INFORMATION
49 int SPflag_;
50 int SPiter_;
51
52 // SECANT INFORMATION
56
57 // TRUNCATED CG INFORMATION
58 Real tol1_;
59 Real tol2_;
60 int maxit_;
61
62 // ALGORITHM SPECIFIC PARAMETERS
63 Real mu0_;
64 Real spexp_;
65 Real alphaMax_;
66
67 mutable int nhess_;
68 unsigned verbosity_;
70
71 bool hasEcon_;
72 Ptr<ReducedLinearConstraint<Real>> rcon_;
73 Ptr<NullSpaceOperator<Real>> ns_;
74
75 using TypeB::Algorithm<Real>::state_;
76 using TypeB::Algorithm<Real>::status_;
77 using TypeB::Algorithm<Real>::proj_;
78
79public:
80 ColemanLiAlgorithm(ParameterList &list, const Ptr<Secant<Real>> &secant = nullPtr);
81
82 using TypeB::Algorithm<Real>::run;
83 void run( Vector<Real> &x,
84 const Vector<Real> &g,
85 Objective<Real> &obj,
87 std::ostream &outStream = std::cout) override;
88
89 void writeHeader( std::ostream& os ) const override;
90
91 void writeName( std::ostream& os ) const override;
92
93 void writeOutput( std::ostream& os, const bool write_header = false ) const override;
94
95private:
97 const Vector<Real> &g,
98 Objective<Real> &obj,
100 std::ostream &outStream = std::cout);
101
102 // Compute the projected step s = P(x + alpha*w) - x
103 // Returns the norm of the projected step s
104 // s -- The projected step upon return
105 // w -- The direction vector w (unchanged)
106 // x -- The anchor vector x (unchanged)
107 // alpha -- The step size (unchanged)
108 // bnd -- The bound constraint
109 Real dgpstep(Vector<Real> &s, const Vector<Real> &w,
110 const Vector<Real> &x, const Real alpha,
111 std::ostream &outStream = std::cout) const;
112
113 // Compute sigma such that ||x+sigma*p||_inv(M) = del. This is called
114 // if dtrpcg detects negative curvature or if the step violates
115 // the trust region bound
116 // xtx -- The dot product <x, inv(M)x> (unchanged)
117 // ptp -- The dot product <p, inv(M)p> (unchanged)
118 // ptx -- The dot product <p, inv(M)x> (unchanged)
119 // del -- Trust region radius (unchanged)
120 Real dtrqsol(const Real xtx, const Real ptp, const Real ptx, const Real del) const;
121
122 // Solve the trust region subproblem: minimize q(w) subject to the
123 // trust region constraint ||w||_inv(M) <= del using the Steihaug-Toint
124 // Conjugate Gradients algorithm
125 // w -- The step vector to be computed
126 // iflag -- Termination flag
127 // iter -- Number of CG iterations
128 // g -- Vector containing gradient (dual space)
129 // x -- Vector containing iterate (primal space)
130 // gdual -- Vector containing primal gradient (primal space)
131 // del -- Trust region radius (unchanged)
132 // model -- Trust region model
133 // bnd -- Bound constraint used to remove active variables
134 // tol -- Residual stopping tolerance (unchanged)
135 // stol -- Preconditioned residual stopping tolerance (unchanged)
136 // p -- Primal working array that stores the CG step
137 // q -- Dual working array that stores the Hessian applied to p
138 // r -- Primal working array that stores the preconditioned residual
139 // t -- Dual working array that stores the residual
140 // pwa1 -- Primal working array that stores the pruned vector in hessVec
141 // pwa2 -- Primal working array that stores the pruned vector in hessVec
142 // dwa -- Dual working array that stores the pruned vector in precond
143 Real dtrpcg(Vector<Real> &w, int &iflag, int &iter,
144 const Vector<Real> &g, const Vector<Real> &x, const Vector<Real> &gdual,
145 const Real del, TrustRegionModel_U<Real> &model, BoundConstraint<Real> &bnd,
146 const Real tol, const Real stol,
148 Vector<Real> &t, Vector<Real> &pwa1, Vector<Real> &pwa2,
149 Vector<Real> &dwa,
150 std::ostream &outStream = std::cout) const;
151
152 void applyC(Vector<Real> &Cv,
153 const Vector<Real> &v,
154 const Vector<Real> &x,
155 const Vector<Real> &g,
157 Vector<Real> &pwa) const;
158
159 void applyHessian(Vector<Real> &hv,
160 const Vector<Real> &v,
161 const Vector<Real> &x,
162 const Vector<Real> &g,
165 Real &tol,
166 Vector<Real> &pwa1,
167 Vector<Real> &pwa2) const;
168
169 void applyPrecond(Vector<Real> &hv,
170 const Vector<Real> &v,
171 const Vector<Real> &x,
172 const Vector<Real> &g,
175 Real &tol,
176 Vector<Real> &dwa,
177 Vector<Real> &pwa) const;
178
179}; // class ROL::TypeB::ColemanLiAlgorithm
180
181} // namespace TypeB
182} // namespace ROL
183
185
186#endif
Provides the interface to apply upper and lower bound constraints.
Provides the interface to evaluate objective functions.
Provides interface for and implements limited-memory secant operators.
Provides the interface to evaluate trust-region model functions.
Provides an interface to run bound constrained optimization algorithms.
Ptr< PolyhedralProjection< Real > > proj_
const Ptr< AlgorithmState< Real > > state_
const Ptr< CombinedStatusTest< Real > > status_
Provides an interface to run the affine-scaling trust-region algorithm of Coleman and Li.
Real dgpstep(Vector< Real > &s, const Vector< Real > &w, const Vector< Real > &x, const Real alpha, std::ostream &outStream=std::cout) const
int nhess_
Number of Hessian applications.
bool useSecantPrecond_
Flag to use secant as a preconditioner (default: false)
int maxit_
Maximum number of CG iterations (default: 20)
void writeOutput(std::ostream &os, const bool write_header=false) const override
Print iterate status.
Ptr< ReducedLinearConstraint< Real > > rcon_
Equality constraint restricted to current active variables.
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout)
bool writeHeader_
Flag to write header at every iteration.
Real dtrpcg(Vector< Real > &w, int &iflag, int &iter, const Vector< Real > &g, const Vector< Real > &x, const Vector< Real > &gdual, const Real del, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, const Real tol, const Real stol, Vector< Real > &p, Vector< Real > &q, Vector< Real > &r, Vector< Real > &t, Vector< Real > &pwa1, Vector< Real > &pwa2, Vector< Real > &dwa, std::ostream &outStream=std::cout) const
Real dtrqsol(const Real xtx, const Real ptp, const Real ptx, const Real del) const
Real alphaMax_
Maximum value of relaxation parameter (default: 0.999)
void applyPrecond(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, Real &tol, Vector< Real > &dwa, Vector< Real > &pwa) const
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...
int SPflag_
Subproblem solver termination flag.
int SPiter_
Subproblem solver iteration count.
TRUtils::ETRFlag TRflag_
Trust-region exit flag.
ESecant esec_
Secant type (default: Limited-Memory BFGS)
bool useSecantHessVec_
Flag to use secant as Hessian (default: false)
unsigned verbosity_
Output level (default: 0)
bool interpRad_
Interpolate the trust-region radius if ratio is negative (default: false)
Real spexp_
Relative tolerance exponent for subproblem solve (default: 1, range: [1,2])
Real gamma0_
Radius decrease rate (negative rho) (default: 0.0625)
Real tol1_
Absolute tolerance for truncated CG (default: 1e-4)
Real eta1_
Radius decrease threshold (default: 0.05)
Real gamma1_
Radius decrease rate (positive rho) (default: 0.25)
Real delMax_
Maximum trust-region radius (default: ROL_INF)
Ptr< NullSpaceOperator< Real > > ns_
Null space projection onto reduced equality constraint Jacobian.
Real eta2_
Radius increase threshold (default: 0.9)
void writeName(std::ostream &os) const override
Print step name.
bool hasEcon_
Flag signifies if equality constraints exist.
Real eta0_
Step acceptance threshold (default: 0.05)
Real TRsafe_
Safeguard size for numerically evaluating ratio (default: 1e2)
Real mu0_
Sufficient decrease parameter (default: 1e-2)
Real gamma2_
Radius increase rate (default: 2.5)
void applyHessian(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, Real &tol, Vector< Real > &pwa1, Vector< Real > &pwa2) const
Real eps_
Safeguard for numerically evaluating ratio.
Ptr< TrustRegionModel_U< Real > > model_
Container for trust-region model.
Real tol2_
Relative tolerance for truncated CG (default: 1e-2)
void applyC(Vector< Real > &Cv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, BoundConstraint< Real > &bnd, Vector< Real > &pwa) const
void writeHeader(std::ostream &os) const override
Print iterate header.
Defines the linear algebra or vector space interface.