ROL
ROL_TruncatedCG.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_TRUNCATEDCG_H
45 #define ROL_TRUNCATEDCG_H
46 
51 #include "ROL_TrustRegion.hpp"
52 #include "ROL_Types.hpp"
53 #include "ROL_HelperFunctions.hpp"
54 
55 namespace ROL {
56 
57 template<class Real>
58 class TruncatedCG : public TrustRegion<Real> {
59 private:
60  Teuchos::RCP<Vector<Real> > primalVector_;
61 
62  Teuchos::RCP<Vector<Real> > s_;
63  Teuchos::RCP<Vector<Real> > g_;
64  Teuchos::RCP<Vector<Real> > v_;
65  Teuchos::RCP<Vector<Real> > p_;
66  Teuchos::RCP<Vector<Real> > Hp_;
67 
68  int maxit_;
69  Real tol1_;
70  Real tol2_;
71 
72  Real pRed_;
73 
74 public:
75 
76  // Constructor
77  TruncatedCG( Teuchos::ParameterList &parlist ) : TrustRegion<Real>(parlist), pRed_(0) {
78  // Unravel Parameter List
79  Real em4(1e-4), em2(1e-2);
80  maxit_ = parlist.sublist("General").sublist("Krylov").get("Iteration Limit",20);
81  tol1_ = parlist.sublist("General").sublist("Krylov").get("Absolute Tolerance",em4);
82  tol2_ = parlist.sublist("General").sublist("Krylov").get("Relative Tolerance",em2);
83  }
84 
85  void initialize( const Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g) {
87 
88  primalVector_ = x.clone();
89 
90  s_ = s.clone();
91  g_ = g.clone();
92  v_ = s.clone();
93  p_ = s.clone();
94  Hp_ = g.clone();
95  }
96 
97  void run( Vector<Real> &s,
98  Real &snorm,
99  int &iflag,
100  int &iter,
101  const Real del,
102  TrustRegionModel<Real> &model ) {
103  Real tol = std::sqrt(ROL_EPSILON<Real>());
104  const Real zero(0), one(1), two(2), half(0.5);
105  // Initialize step
106  s.zero(); s_->zero();
107  snorm = zero;
108  Real snorm2(0), s1norm2(0);
109  // Compute (projected) gradient
110  model.dualTransform(*g_,*model.getGradient());
111  Real gnorm = g_->norm(), normg = gnorm;
112  const Real gtol = std::min(tol1_,tol2_*gnorm);
113  // Preconditioned (projected) gradient vector
114  model.precond(*v_,*g_,s,tol);
115  // Initialize basis vector
116  p_->set(*v_); p_->scale(-one);
117  Real pnorm2 = v_->dot(g_->dual());
118  if ( pnorm2 <= zero ) {
119  iflag = 4;
120  iter = 0;
121  return;
122  }
123  // Initialize scalar storage
124  iter = 0; iflag = 0;
125  Real kappa(0), beta(0), sigma(0), alpha(0), tmp(0), sMp(0);
126  Real gv = v_->dot(g_->dual());
127  pRed_ = zero;
128  // Iterate CG
129  for (iter = 0; iter < maxit_; iter++) {
130  // Apply Hessian to direction p
131  model.hessVec(*Hp_,*p_,s,tol);
132  // Check positivity of Hessian
133  kappa = p_->dot(Hp_->dual());
134  if (kappa <= zero) {
135  sigma = (-sMp+sqrt(sMp*sMp+pnorm2*(del*del-snorm2)))/pnorm2;
136  s.axpy(sigma,*p_);
137  iflag = 2;
138  break;
139  }
140  // Update step
141  alpha = gv/kappa;
142  s_->set(s);
143  s_->axpy(alpha,*p_);
144  s1norm2 = snorm2 + two*alpha*sMp + alpha*alpha*pnorm2;
145  // Check if step exceeds trust region radius
146  if (s1norm2 >= del*del) {
147  sigma = (-sMp+sqrt(sMp*sMp+pnorm2*(del*del-snorm2)))/pnorm2;
148  s.axpy(sigma,*p_);
149  iflag = 3;
150  break;
151  }
152  // Update model predicted reduction
153  pRed_ += half*alpha*gv;
154  // Set step to temporary step and store norm
155  s.set(*s_);
156  snorm2 = s1norm2;
157  // Check for convergence
158  g_->axpy(alpha,*Hp_);
159  normg = g_->norm();
160  if (normg < gtol) {
161  break;
162  }
163  // Preconditioned updated (projected) gradient vector
164  model.precond(*v_,*g_,s,tol);
165  tmp = gv;
166  gv = v_->dot(g_->dual());
167  beta = gv/tmp;
168  // Update basis vector
169  p_->scale(beta);
170  p_->axpy(-one,*v_);
171  sMp = beta*(sMp+alpha*pnorm2);
172  pnorm2 = gv + beta*beta*pnorm2;
173  }
174  // Update model predicted reduction
175  if (iflag > 0) {
176  pRed_ += sigma*(gv-half*sigma*kappa);
177  }
178  // Check iteration count
179  if (iter == maxit_) {
180  iflag = 1;
181  }
182  if (iflag != 1) {
183  iter++;
184  }
185  // Update norm of step and update model predicted reduction
186  model.primalTransform(*s_,s);
187  s.set(*s_);
188  snorm = s.norm();
190  }
191 
192 #if 0
193  void truncatedCG_proj( Vector<Real> &s, Real &snorm, Real &del, int &iflag, int &iter, const Vector<Real> &x,
194  const Vector<Real> &grad, const Real &gnorm, ProjectedObjective<Real> &pObj ) {
195  Real tol = std::sqrt(ROL_EPSILON<Real>());
196 
197  const Real gtol = std::min(tol1_,tol2_*gnorm);
198 
199  // Compute Cauchy Point
200  Real scnorm = 0.0;
201  Teuchos::RCP<Vector<Real> > sc = x.clone();
202  cauchypoint(*sc,scnorm,del,iflag,iter,x,grad,gnorm,pObj);
203  Teuchos::RCP<Vector<Real> > xc = x.clone();
204  xc->set(x);
205  xc->plus(*sc);
206 
207  // Old and New Step Vectors
208  s.set(*sc);
209  snorm = s.norm();
210  Real snorm2 = snorm*snorm;
211  Teuchos::RCP<Vector<Real> > s1 = x.clone();
212  s1->zero();
213  Real s1norm2 = 0.0;
214 
215  // Gradient Vector
216  Teuchos::RCP<Vector<Real> > g = x.clone();
217  g->set(grad);
218  Teuchos::RCP<Vector<Real> > Hs = x.clone();
219  pObj.reducedHessVec(*Hs,s,*xc,x,tol);
220  g->plus(*Hs);
221  Real normg = g->norm();
222 
223  // Preconditioned Gradient Vector
224  Teuchos::RCP<Vector<Real> > v = x.clone();
225  pObj.reducedPrecond(*v,*g,*xc,x,tol);
226 
227  // Basis Vector
228  Teuchos::RCP<Vector<Real> > p = x.clone();
229  p->set(*v);
230  p->scale(-1.0);
231  Real pnorm2 = v->dot(*g);
232 
233  // Hessian Times Basis Vector
234  Teuchos::RCP<Vector<Real> > Hp = x.clone();
235 
236  iter = 0;
237  iflag = 0;
238  Real kappa = 0.0;
239  Real beta = 0.0;
240  Real sigma = 0.0;
241  Real alpha = 0.0;
242  Real tmp = 0.0;
243  Real gv = v->dot(*g);
244  Real sMp = 0.0;
245  pRed_ = 0.0;
246 
247  for (iter = 0; iter < maxit_; iter++) {
248  pObj.reducedHessVec(*Hp,*p,*xc,x,tol);
249 
250  kappa = p->dot(*Hp);
251  if (kappa <= 0) {
252  sigma = (-sMp+sqrt(sMp*sMp+pnorm2*(del*del-snorm2)))/pnorm2;
253  s.axpy(sigma,*p);
254  iflag = 2;
255  break;
256  }
257 
258  alpha = gv/kappa;
259  s1->set(s);
260  s1->axpy(alpha,*p);
261  s1norm2 = snorm2 + 2.0*alpha*sMp + alpha*alpha*pnorm2;
262 
263  if (s1norm2 >= del*del) {
264  sigma = (-sMp+sqrt(sMp*sMp+pnorm2*(del*del-snorm2)))/pnorm2;
265  s.axpy(sigma,*p);
266  iflag = 3;
267  break;
268  }
269 
270  pRed_ += 0.5*alpha*gv;
271 
272  s.set(*s1);
273  snorm2 = s1norm2;
274 
275  g->axpy(alpha,*Hp);
276  normg = g->norm();
277  if (normg < gtol) {
278  break;
279  }
280 
281  pObj.reducedPrecond(*v,*g,*xc,x,tol);
282  tmp = gv;
283  gv = v->dot(*g);
284  beta = gv/tmp;
285 
286  p->scale(beta);
287  p->axpy(-1.0,*v);
288  sMp = beta*(sMp+alpha*pnorm2);
289  pnorm2 = gv + beta*beta*pnorm2;
290  }
291  if (iflag > 0) {
292  pRed_ += sigma*(gv-0.5*sigma*kappa);
293  }
294 
295  if (iter == maxit_) {
296  iflag = 1;
297  }
298  if (iflag != 1) {
299  iter++;
300  }
301 
302  snorm = s.norm();
303  }
304 #endif
305 };
306 
307 }
308 
309 #endif
void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:145
virtual void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
void reducedHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &d, const Vector< Real > &x, Real &tol)
Apply the reduced Hessian to a vector, v. The reduced Hessian first removes elements of v correspondi...
Contains definitions of custom data types in ROL.
Teuchos::RCP< Vector< Real > > primalVector_
Teuchos::RCP< Vector< Real > > p_
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.
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:159
Contains definitions for helper functions in ROL.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:76
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
Teuchos::RCP< Vector< Real > > v_
virtual const Teuchos::RCP< const Vector< Real > > getGradient(void) const
void setPredictedReduction(const Real pRed)
Teuchos::RCP< Vector< Real > > Hp_
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply Hessian approximation to vector.
void run(Vector< Real > &s, Real &snorm, int &iflag, int &iter, const Real del, TrustRegionModel< Real > &model)
Teuchos::RCP< Vector< Real > > g_
Provides interface for truncated CG trust-region subproblem solver.
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:198
virtual Real norm() const =0
Returns where .
void reducedPrecond(Vector< Real > &Mv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &d, const Vector< Real > &x, Real &tol)
Apply the reduced preconditioner to a vector, v. The reduced preconditioner first removes elements of...
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply preconditioner to vector.
Teuchos::RCP< Vector< Real > > s_
TruncatedCG(Teuchos::ParameterList &parlist)
virtual void primalTransform(Vector< Real > &tv, const Vector< Real > &v)