ROL
ROL_SPGTrustRegion_U.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_SPGTRUSTREGION_U_H
45#define ROL_SPGTRUSTREGION_U_H
46
51#include "ROL_TrustRegion_U.hpp"
52#include "ROL_Types.hpp"
53
54#include <deque>
55
56namespace ROL {
57
58template<typename Real>
59class SPGTrustRegion_U : public TrustRegion_U<Real> {
60private:
61 Ptr<Vector<Real>> dwa_, pwa_, pwa1_, gmod_, smin_;
62
65 Real gamma_;
67 int maxit_;
68 Real tol1_;
69 Real tol2_;
70 bool useMin_;
72
73public:
74
75 // Constructor
76 SPGTrustRegion_U( ParameterList &parlist ) {
77 ParameterList &list = parlist.sublist("Step").sublist("Trust Region").sublist("SPG");
78 // Spectral projected gradient parameters
79 lambdaMin_ = list.sublist("Solver").get("Minimum Spectral Step Size", 1e-8);
80 lambdaMax_ = list.sublist("Solver").get("Maximum Spectral Step Size", 1e8);
81 gamma_ = list.sublist("Solver").get("Sufficient Decrease Tolerance", 1e-4);
82 maxSize_ = list.sublist("Solver").get("Maximum Storage Size", 10);
83 maxit_ = list.sublist("Solver").get("Iteration Limit", 25);
84 tol1_ = list.sublist("Solver").get("Absolute Tolerance", 1e-4);
85 tol2_ = list.sublist("Solver").get("Relative Tolerance", 1e-2);
86 useMin_ = list.sublist("Solver").get("Use Smallest Model Iterate", true);
87 useNMSP_ = list.sublist("Solver").get("Use Nonmonotone Search", false);
88 }
89
90 void initialize(const Vector<Real> &x, const Vector<Real> &g) {
91 pwa_ = x.clone();
92 pwa1_ = x.clone();
93 smin_ = x.clone();
94 dwa_ = g.clone();
95 gmod_ = g.clone();
96 }
97
99 Real &snorm,
100 Real &pRed,
101 int &iflag,
102 int &iter,
103 const Real del,
104 TrustRegionModel_U<Real> &model ) {
105 const Real zero(0), half(0.5), one(1), two(2), eps(std::sqrt(ROL_EPSILON<Real>()));
106 Real tol(eps), alpha(1), sHs(0), alphaTmp(1), mmax(0), qmin(0), q(0);
107 Real gnorm(0), ss(0), gs(0);
108 std::deque<Real> mqueue; mqueue.push_back(0);
109 gmod_->set(*model.getGradient());
110
111 // Compute Cauchy point
112 pwa1_->set(gmod_->dual());
113 s.set(*pwa1_); s.scale(-one);
114 model.hessVec(*dwa_,s,s,tol);
115 gs = gmod_->apply(s);
116 sHs = dwa_->apply(s);
117 snorm = std::sqrt(std::abs(gs));
118 alpha = -gs/sHs;
119 if (alpha*snorm >= del || sHs <= zero) alpha = del/snorm;
120 q = alpha*(gs+half*alpha*sHs);
121 gmod_->axpy(alpha,*dwa_);
122 s.scale(alpha);
123
124 if (useNMSP_ && useMin_) { qmin = q; smin_->set(s);}
125
126 // Compute initial projected gradient
127 pwa1_->set(gmod_->dual());
128 pwa_->set(s); pwa_->axpy(-one,*pwa1_);
129 snorm = pwa_->norm();
130 if (snorm > del) pwa_->scale(del/snorm);
131 pwa_->axpy(-one,s);
132 gnorm = pwa_->norm();
133 if (gnorm == zero) {
134 snorm = s.norm();
135 pRed = -q;
136 return;
137 }
138 const Real gtol = std::min(tol1_,tol2_*gnorm);
139
140 // Compute initial step
141 Real lambda = std::max(lambdaMin_,std::min(one/gmod_->norm(),lambdaMax_));
142 pwa_->set(s); pwa_->axpy(-lambda,*pwa1_);
143 snorm = pwa_->norm();
144 if (snorm > del) pwa_->scale(del/snorm);
145 pwa_->axpy(-one,s);
146 gs = gmod_->apply(*pwa_);
147 ss = pwa_->dot(*pwa_);
148
149 for (iter = 0; iter < maxit_; iter++) {
150 // Evaluate model Hessian
151 model.hessVec(*dwa_,*pwa_,s,tol);
152 sHs = dwa_->apply(*pwa_);
153 // Perform line search
154 if (useNMSP_) { // Nonmonotone
155 mmax = *std::max_element(mqueue.begin(),mqueue.end());
156 alphaTmp = (-(one-gamma_)*gs + std::sqrt(std::pow((one-gamma_)*gs,two)-two*sHs*(q-mmax)))/sHs;
157 }
158 else { // Exact
159 alphaTmp = -gs/sHs;
160 }
161 alpha = (sHs > zero ? std::min(one,std::max(zero,alphaTmp)) : one);
162 // Update model quantities
163 q += alpha*(gs+half*alpha*sHs);
164 gmod_->axpy(alpha,*dwa_);
165 s.axpy(alpha,*pwa_);
166 // Update nonmonotone line search information
167 if (useNMSP_) {
168 if (static_cast<int>(mqueue.size())==maxSize_) mqueue.pop_front();
169 mqueue.push_back(q);
170 if (useMin_ && q <= qmin) { qmin = q; smin_->set(s); }
171 }
172 // Compute Projected gradient norm
173 pwa1_->set(gmod_->dual());
174 pwa_->set(s); pwa_->axpy(-one,*pwa1_);
175 snorm = pwa_->norm();
176 if (snorm > del) pwa_->scale(del/snorm);
177 pwa_->axpy(-one,s);
178 gnorm = pwa_->norm();
179 if (gnorm < gtol) break;
180 // Compute new spectral step
181 lambda = (sHs <= eps ? lambdaMax_ : std::max(lambdaMin_,std::min(ss/sHs,lambdaMax_)));
182 pwa_->set(s); pwa_->axpy(-lambda,*pwa1_);
183 snorm = pwa_->norm();
184 if (snorm > del) pwa_->scale(del/snorm);
185 pwa_->axpy(-one,s);
186 gs = gmod_->apply(*pwa_);
187 ss = pwa_->dot(*pwa_);
188 }
189 if (useNMSP_ && useMin_) { q = qmin; s.set(*smin_); }
190 iflag = (iter==maxit_ ? 1 : 0);
191 pRed = -q;
192 snorm = s.norm();
193 }
194};
195
196} // namespace ROL
197
198#endif
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Contains definitions of custom data types in ROL.
Provides interface for truncated CG trust-region subproblem solver.
SPGTrustRegion_U(ParameterList &parlist)
Ptr< Vector< Real > > pwa1_
void solve(Vector< Real > &s, Real &snorm, Real &pRed, int &iflag, int &iter, const Real del, TrustRegionModel_U< Real > &model)
Ptr< Vector< Real > > smin_
Ptr< Vector< Real > > gmod_
void initialize(const Vector< Real > &x, const Vector< Real > &g)
Ptr< Vector< Real > > dwa_
Ptr< Vector< Real > > pwa_
Provides the interface to evaluate trust-region model functions.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Apply Hessian approximation to vector.
virtual const Ptr< const Vector< Real > > getGradient(void) const
Provides interface for and implements trust-region subproblem solvers.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:84
virtual Real norm() const =0
Returns where .
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
virtual void scale(const Real alpha)=0
Compute 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 .
Definition: ROL_Vector.hpp:153