ROL
ROL_TypeU_TrustRegionAlgorithm_Def.hpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Rapid Optimization Library (ROL) Package
5// Copyright (2014) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact lead developers:
38// Drew Kouri (dpkouri@sandia.gov) and
39// Denis Ridzal (dridzal@sandia.gov)
40//
41// ************************************************************************
42// @HEADER
43
44#ifndef ROL_TRUSTREGIONALGORITHM_U_DEF_H
45#define ROL_TRUSTREGIONALGORITHM_U_DEF_H
46
48
49namespace ROL {
50namespace TypeU {
51
52template<typename Real>
54 const Ptr<Secant<Real>> &secant )
55 : Algorithm<Real>(), esec_(SECANT_USERDEFINED) {
56 // Set status test
57 status_->reset();
58 status_->add(makePtr<StatusTest<Real>>(parlist));
59
60 // Trust-Region Parameters
61 ParameterList &slist = parlist.sublist("Step");
62 ParameterList &trlist = slist.sublist("Trust Region");
63 state_->searchSize = trlist.get("Initial Radius", static_cast<Real>(-1));
64 delMax_ = trlist.get("Maximum Radius", ROL_INF<Real>());
65 eta0_ = trlist.get("Step Acceptance Threshold", static_cast<Real>(0.05));
66 eta1_ = trlist.get("Radius Shrinking Threshold", static_cast<Real>(0.05));
67 eta2_ = trlist.get("Radius Growing Threshold", static_cast<Real>(0.9));
68 gamma0_ = trlist.get("Radius Shrinking Rate (Negative rho)", static_cast<Real>(0.0625));
69 gamma1_ = trlist.get("Radius Shrinking Rate (Positive rho)", static_cast<Real>(0.25));
70 gamma2_ = trlist.get("Radius Growing Rate", static_cast<Real>(2.5));
71 TRsafe_ = trlist.get("Safeguard Size", static_cast<Real>(100.0));
72 eps_ = TRsafe_*ROL_EPSILON<Real>();
73 // Inexactness Information
74 ParameterList &glist = parlist.sublist("General");
75 useInexact_.clear();
76 useInexact_.push_back(glist.get("Inexact Objective Function", false));
77 useInexact_.push_back(glist.get("Inexact Gradient", false));
78 useInexact_.push_back(glist.get("Inexact Hessian-Times-A-Vector", false));
79 // Trust-Region Inexactness Parameters
80 ParameterList &ilist = trlist.sublist("Inexact").sublist("Gradient");
81 scale0_ = ilist.get("Tolerance Scaling", static_cast<Real>(0.1));
82 scale1_ = ilist.get("Relative Tolerance", static_cast<Real>(2));
83 // Inexact Function Evaluation Information
84 ParameterList &vlist = trlist.sublist("Inexact").sublist("Value");
85 scale_ = vlist.get("Tolerance Scaling", static_cast<Real>(1.e-1));
86 omega_ = vlist.get("Exponent", static_cast<Real>(0.9));
87 force_ = vlist.get("Forcing Sequence Initial Value", static_cast<Real>(1.0));
88 updateIter_ = vlist.get("Forcing Sequence Update Frequency", static_cast<int>(10));
89 forceFactor_ = vlist.get("Forcing Sequence Reduction Factor", static_cast<Real>(0.1));
90 // Initialize Trust Region Subproblem Solver Object
91 etr_ = StringToETrustRegionU(trlist.get("Subproblem Solver", "Dogleg"));
92 solver_ = TrustRegionUFactory<Real>(parlist);
93 verbosity_ = glist.get("Output Level", 0);
94 // Secant Information
95 useSecantPrecond_ = glist.sublist("Secant").get("Use as Preconditioner", false);
96 useSecantHessVec_ = glist.sublist("Secant").get("Use as Hessian", false);
97 if (secant == nullPtr) {
98 esec_ = StringToESecant(glist.sublist("Secant").get("Type","Limited-Memory BFGS"));
99 }
100 // Initialize trust region model
101 model_ = makePtr<TrustRegionModel_U<Real>>(parlist,secant);
103}
104
105template<typename Real>
107 const Vector<Real> &g,
108 Vector<Real> &Bg,
109 Objective<Real> &obj,
110 std::ostream &outStream) {
111 // Initialize data
113 solver_->initialize(x,g);
114 model_->initialize(x,g);
115 // Update approximate gradient and approximate objective function.
116 Real ftol = static_cast<Real>(0.1)*ROL_OVERFLOW<Real>();
117 obj.update(x,UpdateType::Initial,state_->iter);
118 state_->value = obj.value(x,ftol);
119 state_->nfval++;
120 state_->snorm = ROL_INF<Real>();
121 state_->gnorm = ROL_INF<Real>();
122 computeGradient(x,obj);
123 // Check if inverse Hessian is implemented for dogleg methods
124 model_->validate(obj,x,g,etr_);
125 // Compute initial trust region radius if desired.
126 if ( state_->searchSize <= static_cast<Real>(0) ) {
127 int nfval = 0;
128 state_->searchSize
129 = TRUtils::initialRadius<Real>(nfval,x,*state_->gradientVec,Bg,
130 state_->value,state_->gnorm,obj,*model_,delMax_,
131 outStream,(verbosity_>1));
132 state_->nfval += nfval;
133 }
134}
135
136template<typename Real>
138 Objective<Real> &obj,
139 Real pRed) {
140 const Real one(1);
141 Real tol(std::sqrt(ROL_EPSILON<Real>())), fval(0);
142 if ( useInexact_[0] ) {
143 if ( !(state_->iter%updateIter_) && (state_->iter != 0) ) {
144 force_ *= forceFactor_;
145 }
146 Real eta = static_cast<Real>(0.999)*std::min(eta1_,one-eta2_);
147 tol = scale_*std::pow(eta*std::min(pRed,force_),one/omega_);
148 state_->value = obj.value(*state_->iterateVec,tol);
149 state_->nfval++;
150 }
151 // Evaluate objective function at new iterate
153 fval = obj.value(x,tol);
154 state_->nfval++;
155 return fval;
156}
157
158template<typename Real>
160 Objective<Real> &obj) {
161 if ( useInexact_[1] ) {
162 const Real one(1);
163 Real gtol1 = scale0_*state_->searchSize;
164 Real gtol0 = gtol1 + one;
165 while ( gtol0 > gtol1 ) {
166 obj.gradient(*state_->gradientVec,x,gtol1);
167 state_->gnorm = state_->gradientVec->norm();
168 gtol0 = gtol1;
169 gtol1 = scale0_*std::min(state_->gnorm,state_->searchSize);
170 }
171 }
172 else {
173 Real gtol = std::sqrt(ROL_EPSILON<Real>());
174 obj.gradient(*state_->gradientVec,x,gtol);
175 state_->gnorm = state_->gradientVec->norm();
176 }
177 state_->ngrad++;
178}
179
180template<typename Real>
182 const Vector<Real> &g,
183 Objective<Real> &obj,
184 std::ostream &outStream ) {
185 const Real zero(0);
186 // Initialize trust-region data
187 Real ftrial(0), pRed(0), rho(0);
188 Ptr<Vector<Real>> gvec = g.clone();
189 initialize(x,g,*gvec,obj,outStream);
190
191 // Output
192 if (verbosity_ > 0) writeOutput(outStream,true);
193
194 while (status_->check(*state_)) {
195 // Build trust-region model
196 model_->setData(obj,x,*state_->gradientVec);
197 // Minimize trust-region model over trust-region constraint
198 pRed = zero;
199 SPflag_ = 0; SPiter_ = 0;
200 solver_->solve(*state_->stepVec,state_->snorm,pRed,SPflag_,SPiter_,
201 state_->searchSize,*model_);
202 // Compute trial objective function value
203 x.plus(*state_->stepVec);
204 ftrial = computeValue(x,obj,pRed);
205 // Compute ratio of actual and predicted reduction
206 TRflag_ = TRUtils::SUCCESS;
207 TRUtils::analyzeRatio<Real>(rho,TRflag_,state_->value,ftrial,pRed,eps_,outStream,verbosity_>1);
208 // Update algorithm state
209 state_->iter++;
210 // Accept/reject step and update trust region radius
211 if ((rho < eta0_ && TRflag_ == TRUtils::SUCCESS)
212 || (TRflag_ >= 2)) { // Step Rejected
213 x.set(*state_->iterateVec);
214 obj.update(x,UpdateType::Revert,state_->iter);
215 if (rho < zero && TRflag_ != TRUtils::TRNAN) {
216 // Negative reduction, interpolate to find new trust-region radius
217 state_->searchSize = TRUtils::interpolateRadius<Real>(*state_->gradientVec,*state_->stepVec,
218 state_->snorm,pRed,state_->value,ftrial,state_->searchSize,gamma0_,gamma1_,eta2_,
219 outStream,verbosity_>1);
220 }
221 else { // Shrink trust-region radius
222 state_->searchSize = gamma1_*std::min(state_->snorm,state_->searchSize);
223 }
224 if (useInexact_[1]) computeGradient(x,obj);
225 }
226 else if ((rho >= eta0_ && TRflag_ != TRUtils::NPOSPREDNEG)
227 || (TRflag_ == TRUtils::POSPREDNEG)) { // Step Accepted
228 state_->iterateVec->set(x);
229 state_->value = ftrial;
230 obj.update(x,UpdateType::Accept,state_->iter);
231 // Increase trust-region radius
232 if (rho >= eta2_) state_->searchSize = std::min(gamma2_*state_->searchSize, delMax_);
233 // Compute gradient at new iterate
234 gvec->set(*state_->gradientVec);
235 computeGradient(x,obj);
236 // Update secant information in trust-region model
237 model_->update(x,*state_->stepVec,*gvec,*state_->gradientVec,
238 state_->snorm,state_->iter);
239 }
240 // Update Output
241 if (verbosity_ > 0) writeOutput(outStream,printHeader_);
242 }
243 if (verbosity_ > 0) Algorithm<Real>::writeExitStatus(outStream);
244}
245
246template<typename Real>
247void TrustRegionAlgorithm<Real>::writeHeader( std::ostream& os ) const {
248 std::stringstream hist;
249 if(verbosity_ > 1) {
250 hist << std::string(114,'-') << std::endl;
251 hist << "Trust-Region status output definitions" << std::endl << std::endl;
252 hist << " iter - Number of iterates (steps taken)" << std::endl;
253 hist << " value - Objective function value" << std::endl;
254 hist << " gnorm - Norm of the gradient" << std::endl;
255 hist << " snorm - Norm of the step (update to optimization vector)" << std::endl;
256 hist << " delta - Trust-Region radius" << std::endl;
257 hist << " #fval - Number of times the objective function was evaluated" << std::endl;
258 hist << " #grad - Number of times the gradient was computed" << std::endl;
259 hist << std::endl;
260 hist << " tr_flag - Trust-Region flag" << std::endl;
261 for( int flag = TRUtils::SUCCESS; flag != TRUtils::UNDEFINED; ++flag ) {
262 hist << " " << NumberToString(flag) << " - "
263 << TRUtils::ETRFlagToString(static_cast<TRUtils::ETRFlag>(flag)) << std::endl;
264 }
265 if( etr_ == TRUSTREGION_U_TRUNCATEDCG ) {
266 hist << std::endl;
267 hist << " iterCG - Number of Truncated CG iterations" << std::endl << std::endl;
268 hist << " flagGC - Trust-Region Truncated CG flag" << std::endl;
269 for( int flag = CG_FLAG_SUCCESS; flag != CG_FLAG_UNDEFINED; ++flag ) {
270 hist << " " << NumberToString(flag) << " - "
271 << ECGFlagToString(static_cast<ECGFlag>(flag)) << std::endl;
272 }
273 }
274 else if( etr_ == TRUSTREGION_U_SPG ) {
275 hist << std::endl;
276 hist << " iterCG - Number of spectral projected gradient iterations" << std::endl << std::endl;
277 hist << " flagGC - Trust-Region spectral projected gradient flag" << std::endl;
278 }
279 hist << std::string(114,'-') << std::endl;
280 }
281 hist << " ";
282 hist << std::setw(6) << std::left << "iter";
283 hist << std::setw(15) << std::left << "value";
284 hist << std::setw(15) << std::left << "gnorm";
285 hist << std::setw(15) << std::left << "snorm";
286 hist << std::setw(15) << std::left << "delta";
287 hist << std::setw(10) << std::left << "#fval";
288 hist << std::setw(10) << std::left << "#grad";
289 hist << std::setw(10) << std::left << "tr_flag";
290 if ( etr_ == TRUSTREGION_U_TRUNCATEDCG ) {
291 hist << std::setw(10) << std::left << "iterCG";
292 hist << std::setw(10) << std::left << "flagCG";
293 }
294 else if (etr_ == TRUSTREGION_U_SPG) {
295 hist << std::setw(10) << std::left << "iterSPG";
296 hist << std::setw(10) << std::left << "flagSPG";
297 }
298 hist << std::endl;
299 os << hist.str();
300}
301
302template<typename Real>
303void TrustRegionAlgorithm<Real>::writeName( std::ostream& os ) const {
304 std::stringstream hist;
305 hist << std::endl << ETrustRegionUToString(etr_) << " Trust-Region Solver";
306 if ( useSecantPrecond_ || useSecantHessVec_ ) {
307 if ( useSecantPrecond_ && !useSecantHessVec_ ) {
308 hist << " with " << ESecantToString(esec_) << " Preconditioning" << std::endl;
309 }
310 else if ( !useSecantPrecond_ && useSecantHessVec_ ) {
311 hist << " with " << ESecantToString(esec_) << " Hessian Approximation" << std::endl;
312 }
313 else {
314 hist << " with " << ESecantToString(esec_) << " Preconditioning and Hessian Approximation" << std::endl;
315 }
316 }
317 else {
318 hist << std::endl;
319 }
320 os << hist.str();
321}
322
323template<typename Real>
324void TrustRegionAlgorithm<Real>::writeOutput(std::ostream& os, bool print_header) const {
325 std::stringstream hist;
326 hist << std::scientific << std::setprecision(6);
327 if ( state_->iter == 0 ) {
328 writeName(os);
329 }
330 if ( print_header ) {
331 writeHeader(os);
332 }
333 if ( state_->iter == 0 ) {
334 hist << " ";
335 hist << std::setw(6) << std::left << state_->iter;
336 hist << std::setw(15) << std::left << state_->value;
337 hist << std::setw(15) << std::left << state_->gnorm;
338 hist << std::setw(15) << std::left << "---";
339 hist << std::setw(15) << std::left << state_->searchSize;
340 hist << std::setw(10) << std::left << state_->nfval;
341 hist << std::setw(10) << std::left << state_->ngrad;
342 hist << std::setw(10) << std::left << "---";
343 if ( etr_ == TRUSTREGION_U_TRUNCATEDCG || etr_ == TRUSTREGION_U_SPG ) {
344 hist << std::setw(10) << std::left << "---";
345 hist << std::setw(10) << std::left << "---";
346 }
347 hist << std::endl;
348 }
349 else {
350 hist << " ";
351 hist << std::setw(6) << std::left << state_->iter;
352 hist << std::setw(15) << std::left << state_->value;
353 hist << std::setw(15) << std::left << state_->gnorm;
354 hist << std::setw(15) << std::left << state_->snorm;
355 hist << std::setw(15) << std::left << state_->searchSize;
356 hist << std::setw(10) << std::left << state_->nfval;
357 hist << std::setw(10) << std::left << state_->ngrad;
358 hist << std::setw(10) << std::left << TRflag_;
359 if ( etr_ == TRUSTREGION_U_TRUNCATEDCG || etr_ == TRUSTREGION_U_SPG ) {
360 hist << std::setw(10) << std::left << SPiter_;
361 hist << std::setw(10) << std::left << SPflag_;
362 }
363 hist << std::endl;
364 }
365 os << hist.str();
366}
367} // namespace TypeU
368} // namespace ROL
369
370#endif
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Contains definitions of enums for trust region algorithms.
Provides the interface to evaluate objective functions.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
virtual void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
Provides interface for and implements limited-memory secant operators.
Definition: ROL_Secant.hpp:79
Provides an interface to check status of optimization algorithms.
Provides an interface to run unconstrained optimization algorithms.
const Ptr< CombinedStatusTest< Real > > status_
void initialize(const Vector< Real > &x, const Vector< Real > &g)
virtual void writeExitStatus(std::ostream &os) const
const Ptr< AlgorithmState< Real > > state_
Real scale1_
Scale for inexact gradient computation.
bool printHeader_
Print header at every iteration.
std::vector< bool > useInexact_
Flags for inexact (0) objective function, (1) gradient, (2) Hessian.
Real gamma0_
Radius decrease rate (negative rho).
Real TRsafe_
Safeguard size for numerically evaluating ratio.
void writeHeader(std::ostream &os) const override
Print iterate header.
Real scale0_
Scale for inexact gradient computation.
void writeOutput(std::ostream &os, bool print_header=false) const override
Print iterate status.
Real eps_
Safeguard for numerically evaluating ratio.
int verbosity_
Print additional information to screen if > 0.
void writeName(std::ostream &os) const override
Print step name.
Real computeValue(const Vector< Real > &x, Objective< Real > &obj, Real pRed)
Real gamma1_
Radius decrease rate (positive rho).
Real delMax_
Maximum trust-region radius.
Ptr< TrustRegionModel_U< Real > > model_
Container for trust-region model.
TrustRegionAlgorithm(ParameterList &parlist, const Ptr< Secant< Real > > &secant=nullPtr)
void initialize(const Vector< Real > &x, const Vector< Real > &g, Vector< Real > &Bg, Objective< Real > &obj, std::ostream &outStream=std::cout)
void run(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, std::ostream &outStream=std::cout) override
Run algorithm on unconstrained problems (Type-U). This general interface supports the use of dual opt...
ETrustRegionU etr_
Trust-region subproblem solver type.
void computeGradient(const Vector< Real > &x, Objective< Real > &obj)
Compute gradient to iteratively satisfy inexactness condition.
Ptr< TrustRegion_U< Real > > solver_
Container for trust-region solver object.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:84
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
virtual void plus(const Vector &x)=0
Compute , where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
std::string ETRFlagToString(ETRFlag trf)
std::string NumberToString(T Number)
Definition: ROL_Types.hpp:81
ESecant StringToESecant(std::string s)
Definition: ROL_Types.hpp:541
@ SECANT_USERDEFINED
Definition: ROL_Types.hpp:489
ETrustRegionU StringToETrustRegionU(std::string s)
@ CG_FLAG_UNDEFINED
Definition: ROL_Types.hpp:825
@ CG_FLAG_SUCCESS
Definition: ROL_Types.hpp:820
std::string ESecantToString(ESecant tr)
Definition: ROL_Types.hpp:493
std::string ECGFlagToString(ECGFlag cgf)
Definition: ROL_Types.hpp:829
std::string ETrustRegionUToString(ETrustRegionU tr)