ROL
ROL_TypeP_InexactNewtonAlgorithm_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_TYPEP_QUASINEWTONALGORITHM_DEF_HPP
45#define ROL_TYPEP_QUASINEWTONALGORITHM_DEF_HPP
46
51
52namespace ROL {
53namespace TypeP {
54
55template<typename Real>
57 : list_(list) {
58 // Set status test
59 status_->reset();
60 status_->add(makePtr<StatusTest<Real>>(list));
61
62 // Parse parameter list
63 ParameterList &lslist = list.sublist("Step").sublist("Line Search");
64 t0_ = list.sublist("Status Test").get("Gradient Scale" , 1.0);
65 initProx_ = lslist.get("Apply Prox to Initial Guess", false);
66 maxit_ = lslist.get("Function Evaluation Limit", 20);
67 c1_ = lslist.get("Sufficient Decrease Tolerance", 1e-4);
68 rhodec_ = lslist.sublist("Line-Search Method").get("Backtracking Rate", 0.5);
69 sigma1_ = lslist.sublist("Inexact Newton").get("Lower Step Size Safeguard", 0.1);
70 sigma2_ = lslist.sublist("Inexact Newton").get("Upper Step Size Safeguard", 0.9);
71 algoName_ = lslist.sublist("Inexact Newton").get("Subproblem Solver","Spectral Gradient");
72 int sp_maxit = lslist.sublist("Inexact Newton").get("Subproblem Iteration Limit", 1000);
73 sp_tol1_ = lslist.sublist("Inexact Newton").get("Subproblem Absolute Tolerance", 1e-4);
74 sp_tol2_ = lslist.sublist("Inexact Newton").get("Subproblem Relative Tolerance", 1e-2);
75 sp_exp_ = lslist.sublist("Inexact Newton").get("Subproblem Tolerance Exponent", 1.0);
76 Real opt_tol = lslist.sublist("Status Test").get("Gradient Tolerance", 1e-8);
77 sp_tol_min_ = static_cast<Real>(1e-4)*opt_tol;
78 verbosity_ = list.sublist("General").get("Output Level", 0);
80
81 list_.sublist("Status Test").set("Iteration Limit", sp_maxit);
82 list_.sublist("General").set("Output Level", verbosity_>0 ? verbosity_-1 : 0);
83}
84
85
86template<typename Real>
88 const Vector<Real> &g,
89 Objective<Real> &sobj,
90 Objective<Real> &nobj,
91 Vector<Real> &dg,
92 Vector<Real> &px,
93 std::ostream &outStream) {
94 const Real one(1);
95 Real tol(std::sqrt(ROL_EPSILON<Real>()));
96 // Initialize data
98 // Update approximate gradient and approximate objective function.
99 Real ftol = std::sqrt(ROL_EPSILON<Real>());
100 if (initProx_) {
101 state_->iterateVec->set(x);
102 nobj.prox(x,*state_->iterateVec,one,tol); state_->nprox++;
103 }
104 sobj.update(x,UpdateType::Initial,state_->iter);
105 nobj.update(x,UpdateType::Initial,state_->iter);
106 state_->svalue = sobj.value(x,ftol); state_->nsval++;
107 state_->nvalue = nobj.value(x,ftol); state_->nnval++;
108 state_->value = state_->svalue + state_->nvalue;
109 sobj.gradient(*state_->gradientVec,x,ftol); state_->ngrad++;
110 dg.set(state_->gradientVec->dual());
111 pgstep(*state_->iterateVec,px,nobj,x,dg,t0_,tol);
112 state_->gnorm = px.norm() / t0_;
113 state_->snorm = ROL_INF<Real>();
114 nhess_ = 0;
115}
116
117template<typename Real>
119 const Vector<Real> &g,
120 Objective<Real> &sobj,
121 Objective<Real> &nobj,
122 std::ostream &outStream ) {
123 const Real half(0.5), one(1), eps(ROL_EPSILON<Real>());
124 // Initialize trust-region data
125 Ptr<Vector<Real>> s = x.clone(), gp = x.clone(), xs = x.clone(), px = x.clone();
126 initialize(x,g,sobj,nobj,*gp,*px,outStream);
127 Real strial(0), ntrial(0), ftrial(0), gs(0), Qk(0), rhoTmp(0);
128 Real tol(std::sqrt(ROL_EPSILON<Real>())), gtol(1);
129
130 Ptr<TypeP::Algorithm<Real>> algo;
131 Ptr<NewtonObj> qobj = makePtr<NewtonObj>(makePtrFromRef(sobj),x,g);
132
133 // Output
134 if (verbosity_ > 0) writeOutput(outStream,true);
135
136 // Compute steepest descent step
137 xs->set(*state_->iterateVec);
138 state_->iterateVec->set(x);
139 while (status_->check(*state_)) {
140 qobj->setData(x,*state_->gradientVec);
141 // Compute step
142 gtol = std::max(sp_tol_min_,std::min(sp_tol1_,sp_tol2_*std::pow(state_->gnorm,sp_exp_)));
143 list_.sublist("Status Test").set("Gradient Tolerance",gtol);
144 if (algoName_ == "Line Search") algo = makePtr<TypeP::ProxGradientAlgorithm<Real>>(list_);
145 else if (algoName_ == "iPiano") algo = makePtr<TypeP::iPianoAlgorithm<Real>>(list_);
146 else if (algoName_ == "Trust Region") algo = makePtr<TypeP::TrustRegionAlgorithm<Real>>(list_);
147 else algo = makePtr<TypeP::SpectralGradientAlgorithm<Real>>(list_);
148 algo->run(*xs,*qobj,nobj,outStream);
149 s->set(*xs); s->axpy(-one,x);
150 spgIter_ = algo->getState()->iter;
151 nhess_ += qobj->numHessVec();
152 state_->nprox += staticPtrCast<const TypeP::AlgorithmState<Real>>(algo->getState())->nprox;
153
154 // Perform backtracking line search
155 state_->searchSize = one;
156 x.set(*state_->iterateVec);
157 x.axpy(state_->searchSize,*s);
160 strial = sobj.value(x,tol);
161 ntrial = nobj.value(x,tol);
162 ftrial = strial + ntrial;
163 ls_nfval_ = 1;
164 gs = state_->gradientVec->apply(*s);
165 Qk = gs + ntrial - state_->nvalue;
166 if (verbosity_ > 1) {
167 outStream << " In TypeP::InexactNewtonAlgorithm: Line Search" << std::endl;
168 outStream << " Step size: " << state_->searchSize << std::endl;
169 outStream << " Trial objective value: " << ftrial << std::endl;
170 outStream << " Computed reduction: " << state_->value-ftrial << std::endl;
171 outStream << " Dot product of gradient and step: " << gs << std::endl;
172 outStream << " Sufficient decrease bound: " << -Qk*c1_ << std::endl;
173 outStream << " Number of function evaluations: " << ls_nfval_ << std::endl;
174 }
175 if (Qk > -eps) {
176 s->set(*px);
177 x.set(*state_->iterateVec);
178 x.axpy(state_->searchSize,*s);
181 strial = sobj.value(x,tol);
182 ntrial = nobj.value(x,tol);
183 ftrial = strial + ntrial;
184 ls_nfval_++;
185 gs = state_->gradientVec->apply(*s);
186 Qk = gs + ntrial - state_->nvalue;
187 }
188 while ( ftrial > state_->value + c1_*Qk && ls_nfval_ < maxit_ ) {
189 rhoTmp = -half * Qk / (strial-state_->svalue-state_->searchSize*gs);
190 state_->searchSize = ((sigma1_ <= rhoTmp && rhoTmp <= sigma2_) ? rhoTmp : rhodec_) * state_->searchSize;
191 x.set(*state_->iterateVec);
192 x.axpy(state_->searchSize,*s);
195 strial = sobj.value(x,tol);
196 ntrial = nobj.value(x,tol);
197 ftrial = strial + ntrial;
198 Qk = state_->searchSize * gs + ntrial - state_->nvalue;
199 ls_nfval_++;
200 if (verbosity_ > 1) {
201 outStream << std::endl;
202 outStream << " Step size: " << state_->searchSize << std::endl;
203 outStream << " Trial objective value: " << ftrial << std::endl;
204 outStream << " Computed reduction: " << state_->value-ftrial << std::endl;
205 outStream << " Dot product of gradient and step: " << gs << std::endl;
206 outStream << " Sufficient decrease bound: " << -Qk*c1_ << std::endl;
207 outStream << " Number of function evaluations: " << ls_nfval_ << std::endl;
208 }
209 }
210 state_->nsval += ls_nfval_;
211 state_->nnval += ls_nfval_;
212
213 // Compute norm of step
214 state_->stepVec->set(*s);
215 state_->stepVec->scale(state_->searchSize);
216 state_->snorm = state_->stepVec->norm();
217
218 // Update iterate
219 state_->iterateVec->set(x);
220
221 // Compute new value and gradient
222 state_->iter++;
223 state_->value = ftrial;
224 state_->svalue = strial;
225 state_->nvalue = ntrial;
226 sobj.update(x,UpdateType::Accept,state_->iter);
227 nobj.update(x,UpdateType::Accept,state_->iter);
228 sobj.gradient(*state_->gradientVec,x,tol); state_->ngrad++;
229 gp->set(state_->gradientVec->dual());
230
231 // Compute projected gradient norm
232 pgstep(*xs,*px,nobj,x,*gp,t0_,tol);
233 state_->gnorm = s->norm() / t0_;
234
235 // Update Output
236 if (verbosity_ > 0) writeOutput(outStream,writeHeader_);
237 }
239}
240
241template<typename Real>
242void InexactNewtonAlgorithm<Real>::writeHeader( std::ostream& os ) const {
243 std::stringstream hist;
244 if (verbosity_ > 1) {
245 hist << std::string(114,'-') << std::endl;
246 hist << "Line-Search Inexact Proximal Newton";
247 hist << " status output definitions" << std::endl << std::endl;
248 hist << " iter - Number of iterates (steps taken)" << std::endl;
249 hist << " value - Objective function value" << std::endl;
250 hist << " gnorm - Norm of the gradient" << std::endl;
251 hist << " snorm - Norm of the step (update to optimization vector)" << std::endl;
252 hist << " alpha - Line search step length" << std::endl;
253 hist << " #sval - Cumulative number of times the smooth objective function was evaluated" << std::endl;
254 hist << " #nval - Cumulative number of times the nonsmooth objective function was evaluated" << std::endl;
255 hist << " #grad - Cumulative number of times the gradient was computed" << std::endl;
256 hist << " #hess - Cumulative number of times the Hessian was applied" << std::endl;
257 hist << " #prox - Cumulative number of times the projection was computed" << std::endl;
258 hist << " ls_#fval - Number of the times the objective function was evaluated during the line search" << std::endl;
259 hist << " sp_iter - Number iterations to compute quasi-Newton step" << std::endl;
260 hist << std::string(114,'-') << std::endl;
261 }
262
263 hist << " ";
264 hist << std::setw(6) << std::left << "iter";
265 hist << std::setw(15) << std::left << "value";
266 hist << std::setw(15) << std::left << "gnorm";
267 hist << std::setw(15) << std::left << "snorm";
268 hist << std::setw(15) << std::left << "alpha";
269 hist << std::setw(10) << std::left << "#sval";
270 hist << std::setw(10) << std::left << "#nval";
271 hist << std::setw(10) << std::left << "#grad";
272 hist << std::setw(10) << std::left << "#hess";
273 hist << std::setw(10) << std::left << "#prox";
274 hist << std::setw(10) << std::left << "#ls_fval";
275 hist << std::setw(10) << std::left << "sp_iter";
276 hist << std::endl;
277 os << hist.str();
278}
279
280template<typename Real>
281void InexactNewtonAlgorithm<Real>::writeName( std::ostream& os ) const {
282 std::stringstream hist;
283 hist << std::endl << "Line-Search Inexact Proximal Newton (Type P)" << std::endl;
284 os << hist.str();
285}
286
287template<typename Real>
288void InexactNewtonAlgorithm<Real>::writeOutput( std::ostream& os, bool write_header ) const {
289 std::stringstream hist;
290 hist << std::scientific << std::setprecision(6);
291 if ( state_->iter == 0 ) writeName(os);
292 if ( write_header ) writeHeader(os);
293 if ( state_->iter == 0 ) {
294 hist << " ";
295 hist << std::setw(6) << std::left << state_->iter;
296 hist << std::setw(15) << std::left << state_->value;
297 hist << std::setw(15) << std::left << state_->gnorm;
298 hist << std::setw(15) << std::left << "---";
299 hist << std::setw(15) << std::left << "---";
300 hist << std::setw(10) << std::left << state_->nsval;
301 hist << std::setw(10) << std::left << state_->nnval;
302 hist << std::setw(10) << std::left << state_->ngrad;
303 hist << std::setw(10) << std::left << nhess_;
304 hist << std::setw(10) << std::left << state_->nprox;
305 hist << std::setw(10) << std::left << "---";
306 hist << std::setw(10) << std::left << "---";
307 hist << std::endl;
308 }
309 else {
310 hist << " ";
311 hist << std::setw(6) << std::left << state_->iter;
312 hist << std::setw(15) << std::left << state_->value;
313 hist << std::setw(15) << std::left << state_->gnorm;
314 hist << std::setw(15) << std::left << state_->snorm;
315 hist << std::setw(15) << std::left << state_->searchSize;
316 hist << std::setw(10) << std::left << state_->nsval;
317 hist << std::setw(10) << std::left << state_->nnval;
318 hist << std::setw(10) << std::left << state_->ngrad;
319 hist << std::setw(10) << std::left << nhess_;
320 hist << std::setw(10) << std::left << state_->nprox;
321 hist << std::setw(10) << std::left << ls_nfval_;
322 hist << std::setw(10) << std::left << spgIter_;
323 hist << std::endl;
324 }
325 os << hist.str();
326}
327
328} // namespace TypeP
329} // namespace ROL
330
331#endif
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.
Provides the interface to evaluate objective functions.
virtual void prox(Vector< Real > &Pv, const Vector< Real > &v, Real t, Real &tol)
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 an interface to check status of optimization algorithms.
void pgstep(Vector< Real > &pgiter, Vector< Real > &pgstep, Objective< Real > &nobj, const Vector< Real > &x, const Vector< Real > &dg, Real t, Real &tol) const
const Ptr< AlgorithmState< Real > > state_
virtual void writeExitStatus(std::ostream &os) const
const Ptr< CombinedStatusTest< Real > > status_
void initialize(const Vector< Real > &x, const Vector< Real > &g)
Real c1_
Sufficient Decrease Parameter (default: 1e-4).
void run(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &sobj, Objective< Real > &nobj, std::ostream &outStream=std::cout) override
Run algorithm on unconstrained problems (Type-U). This general interface supports the use of dual opt...
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &sobj, Objective< Real > &nobj, Vector< Real > &dg, Vector< Real > &px, std::ostream &outStream=std::cout)
Real sigma2_
Upper safeguard for quadratic line search (default: 0.9).
Real sigma1_
Lower safeguard for quadratic line search (default: 0.1).
void writeName(std::ostream &os) const override
Print step name.
void writeOutput(std::ostream &os, bool write_header=false) const override
Print iterate status.
int maxit_
Maximum number of line search steps (default: 20).
Real rhodec_
Backtracking rate (default: 0.5).
void writeHeader(std::ostream &os) const override
Print iterate header.
Defines the linear algebra or vector space interface.
virtual Real norm() const =0
Returns where .
virtual void set(const Vector &x)
Set where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
Definition ROL_Types.hpp:91
Real ROL_INF(void)