ROL
ROL_CauchyPoint.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_CAUCHYPOINT_H
45 #define ROL_CAUCHYPOINT_H
46 
51 #include "ROL_TrustRegion.hpp"
52 #include "ROL_Vector.hpp"
53 #include "ROL_Types.hpp"
54 #include "ROL_HelperFunctions.hpp"
55 #include "Teuchos_ParameterList.hpp"
56 
57 namespace ROL {
58 
59 template<class Real>
60 class CauchyPoint : public TrustRegion<Real> {
61 private:
62 
63  Teuchos::RCP<Vector<Real> > g_;
64  Teuchos::RCP<Vector<Real> > p_;
65  Teuchos::RCP<Vector<Real> > Hp_;
66 
67  Real pRed_;
68  Real eps_;
69  Real alpha_;
70 
71  bool useCGTCP_;
72 
73 public:
74 
75  // Constructor
76  CauchyPoint( Teuchos::ParameterList &parlist )
77  : TrustRegion<Real>(parlist), pRed_(0), alpha_(-1), useCGTCP_(false) {
78  // Unravel Parameter List
79  Real oe2(100);
80  Real TRsafe = parlist.sublist("Step").sublist("Trust Region").get("Safeguard Size",oe2);
81  eps_ = TRsafe*ROL_EPSILON<Real>();
82  }
83 
84  void initialize( const Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g) {
86  Hp_ = g.clone();
87  p_ = s.clone();
88 // if ( useCGTCP_ ) {
89 // g_ = g.clone();
90 // p_ = s.clone();
91 // }
92  }
93 
94  void run( Vector<Real> &s,
95  Real &snorm,
96  int &iflag,
97  int &iter,
98  const Real del,
99  TrustRegionModel<Real> &model) {
100  //if ( pObj.isConActivated() ) {
101  // if ( useCGTCP_ ) {
102  // cauchypoint_CGT( s, snorm, iflag, iter, del, model );
103  // }
104  // else {
105  // cauchypoint_M( s, snorm, iflag, iter, del, model );
106  // }
107  //}
108  //else {
109  // cauchypoint_unc( s, snorm, iflag, iter, del, model );
110  //}
111  cauchypoint_unc( s, snorm, iflag, iter, del, model );
113  }
114 
115 private:
117  Real &snorm,
118  int &iflag,
119  int &iter,
120  const Real del,
121  TrustRegionModel<Real> &model) {
122  Real tol = std::sqrt(ROL_EPSILON<Real>());
123  // Set step to (projected) gradient
124  model.dualTransform(*Hp_,*model.getGradient());
125  s.set(Hp_->dual());
126  // Apply (reduced) Hessian to (projected) gradient
127  model.hessVec(*Hp_,s,s,tol);
128  Real gBg = Hp_->dot(s.dual());
129  Real gnorm = s.dual().norm();
130  Real gg = gnorm*gnorm;
131  Real alpha = del/gnorm;
132  if ( gBg > ROL_EPSILON<Real>() ) {
133  alpha = std::min(gg/gBg, del/gnorm);
134  }
135 
136  s.scale(-alpha);
137  model.primalTransform(*p_,s);
138  s.set(*p_);
139  snorm = s.norm(); //alpha*gnorm;
140  iflag = 0;
141  iter = 0;
142  pRed_ = alpha*(gg - static_cast<Real>(0.5)*alpha*gBg);
143  }
144 
145 // void cauchypoint_M( Vector<Real> &s,
146 // Real &snorm,
147 // int &iflag,
148 // int &iter,
149 // const Real del,
150 // const Vector<Real> &x,
151 // TrustRegionModel<Real> &model,
152 // BoundConstraint<Real> &bnd) {
153 // Real tol = std::sqrt(ROL_EPSILON<Real>()),
154 // const Real zero(0), half(0.5), oe4(1.e4), two(2);
155 // // Parameters
156 // Real mu0(1.e-2), mu1(1), beta1(0), beta2(0);
157 // bool decr = true;
158 // bool stat = true;
159 // // Initial step length
160 // Real alpha = (alpha_ > zero ? alpha_ : one);
161 // Real alpha0 = alpha;
162 // Real alphamax = oe4*alpha;
163 // // Set step to (projected) gradient
164 // s.zero();
165 // model.gradient(*Hp_,s,tol);
166 // s.set(Hp_->dual());
167 // // Initial model value
168 // s.scale(-alpha);
169 // bnd.computeProjectedStep(s,x);
170 // snorm = s.norm();
171 // Real val = model.value(s,tol);
172 // Real val0 = val;
173 //
174 // // Determine whether to increase or decrease alpha
175 // if ( val > mu0 * gs || snorm > mu1 * del ) {
176 // beta1 = half;
177 // beta2 = half;
178 // decr = true;
179 // }
180 // else {
181 // beta1 = two;
182 // beta2 = two;
183 // decr = false;
184 // }
185 //
186 // while ( stat ) {
187 // // Update step length
188 // alpha0 = alpha;
189 // val0 = val;
190 // alpha *= half*(beta1+beta2);
191 //
192 // // Update model value
193 // s.set(grad.dual());
194 // s.scale(-alpha);
195 // pObj.computeProjectedStep(s,x);
196 // snorm = s.norm();
197 // pObj.hessVec(*Hp_,s,x,tol);
198 // gs = s.dot(grad.dual());
199 // val = gs + half*s.dot(Hp_->dual());
200 //
201 // // Update termination criterion
202 // if ( decr ) {
203 // stat = ( val > mu0 * gs || snorm > mu1 * del );
204 // if ( std::abs(val) < eps_ && std::abs(mu0 *gs) < eps_ ) {
205 // stat = (snorm > mu1 * del);
206 // }
207 // }
208 // else {
209 // stat = !( val > mu0 * gs || snorm > mu1 * del );
210 // if ( std::abs(val) < eps_ && std::abs(mu0 *gs) < eps_ ) {
211 // stat = !(snorm > mu1 * del);
212 // }
213 // if ( alpha > alphamax ) {
214 // stat = false;
215 // }
216 // }
217 // }
218 // // Reset to last 'successful' step
219 // val = val0;
220 // alpha = alpha0;
221 // s.set(grad.dual());
222 // s.scale(-alpha);
223 // pObj.computeProjectedStep(s,x);
224 // snorm = s.norm();
225 //
226 // alpha_ = alpha;
227 // pRed_ = -val;
228 // }
229 //
230 // void cauchypoint_CGT( Vector<Real> &s, Real &snorm, Real &del, int &iflag, int &iter, const Vector<Real> &x,
231 // const Vector<Real> &grad, const Real &gnorm, ProjectedObjective<Real> &pObj ) {
232 // Real tol = std::sqrt(ROL_EPSILON<Real>()), one(1), half(0.5), two(2);
233 // bool tmax_flag = true;
234 // int maxit = 20;
235 // Real t = del/gnorm;
236 // Real tmax(1.e10), tmin(0), gs(0), pgnorm(0);
237 // Real c1(0.25), c2(0.75), c3(0.9), c4(0.25);
238 // for ( int i = 0; i < maxit; i++ ) {
239 // // Compute p = x + s = P(x - t*g)
240 // p_->set(x);
241 // p_->axpy(-t,grad.dual());
242 // pObj.project(*p_);
243 // // Compute s = p - x = P(x - t*g) - x
244 // s.set(*p_);
245 // s.axpy(-one,x);
246 // snorm = s.norm();
247 // // Evaluate Model
248 // pObj.hessVec(*Hp_,s,x,tol);
249 // gs = s.dot(grad.dual());
250 // pRed_ = -gs - half*s.dot(Hp_->dual());
251 //
252 // // Check Stopping Conditions
253 // g_->set(grad);
254 // pObj.pruneActive(*g_,grad,*p_); // Project gradient onto tangent cone at p
255 // pgnorm = g_->norm();
256 // if ( snorm > del || pRed_ < -c2*gs ) {
257 // tmax = t;
258 // tmax_flag = false;
259 // }
260 // else if ( snorm < c3*del && pRed_ > -c1*gs && pgnorm > c4*std::abs(gs)/del ) {
261 // tmin = t;
262 // }
263 // else {
264 // break;
265 // }
266 //
267 // // Update t
268 // if ( tmax_flag ) {
269 // t *= two;
270 // }
271 // else {
272 // t = half*(tmax + tmin);
273 // }
274 // }
275 // }
276 };
277 
278 }
279 
280 #endif
virtual void scale(const Real alpha)=0
Compute where .
CauchyPoint(Teuchos::ParameterList &parlist)
virtual void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
Teuchos::RCP< Vector< Real > > p_
Contains definitions of custom data types in ROL.
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
Provides interface for and implements trust-region subproblem solvers.
Provides the interface to evaluate trust-region model functions.
Contains definitions for helper functions in ROL.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:74
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: ROL_Vector.hpp:213
virtual const Teuchos::RCP< const Vector< Real > > getGradient(void) const
void setPredictedReduction(const Real pRed)
Teuchos::RCP< Vector< Real > > Hp_
void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
Teuchos::RCP< Vector< Real > > g_
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply Hessian approximation to vector.
Provides interface for the Cauchy point trust-region subproblem solver.
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:196
virtual Real norm() const =0
Returns where .
void cauchypoint_unc(Vector< Real > &s, Real &snorm, int &iflag, int &iter, const Real del, TrustRegionModel< Real > &model)
void run(Vector< Real > &s, Real &snorm, int &iflag, int &iter, const Real del, TrustRegionModel< Real > &model)
virtual void primalTransform(Vector< Real > &tv, const Vector< Real > &v)