ROL
ROL_ObjectiveDef.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_OBJECTIVE_DEF_H
45 #define ROL_OBJECTIVE_DEF_H
46 
51 namespace ROL {
52 
53 template <class Real>
54 Real Objective<Real>::dirDeriv( const Vector<Real> &x, const Vector<Real> &d, Real &tol) {
55  Real ftol = std::sqrt(ROL_EPSILON);
56 
57  Teuchos::RCP<Vector<Real> > xd = d.clone();
58  xd->set(x);
59  xd->axpy(tol, d);
60  return (this->value(*xd,ftol) - this->value(x,ftol)) / tol;
61 }
62 
63 template <class Real>
64 void Objective<Real>::gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
65  g.zero();
66  Real deriv = 0.0;
67  Real h = 0.0;
68  for (int i = 0; i < g.dimension(); i++) {
69  h = x.dot(*x.basis(i))*tol;
70  deriv = this->dirDeriv(x,*x.basis(i),h);
71  g.axpy(deriv,*g.basis(i));
72  }
73 }
74 
75 template <class Real>
76 void Objective<Real>::hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
77  Real gtol = std::sqrt(ROL_EPSILON);
78 
79  // Get Step Length
80  Real h = std::max(1.0,x.norm()/v.norm())*tol;
81  //Real h = 2.0/(v.norm()*v.norm())*tol;
82 
83  // Compute Gradient at x
84  Teuchos::RCP<Vector<Real> > g = hv.clone();
85  this->gradient(*g,x,gtol);
86 
87  // Compute New Step x + h*v
88  Teuchos::RCP<Vector<Real> > xnew = x.clone();
89  xnew->set(x);
90  xnew->axpy(h,v);
91  this->update(*xnew);
92 
93  // Compute Gradient at x + h*v
94  hv.zero();
95  this->gradient(hv,*xnew,gtol);
96 
97  // Compute Newton Quotient
98  hv.axpy(-1.0,*g);
99  hv.scale(1.0/h);
100 }
101 
102 
103 
104 template <class Real>
105 std::vector<std::vector<Real> > Objective<Real>::checkGradient( const Vector<Real> &x,
106  const Vector<Real> &g,
107  const Vector<Real> &d,
108  const bool printToStream,
109  std::ostream & outStream,
110  const int numSteps,
111  const int order ) {
112 
113  std::vector<Real> steps(numSteps);
114  for(int i=0;i<numSteps;++i) {
115  steps[i] = pow(10,-i);
116  }
117 
118  return checkGradient(x,g,d,steps,printToStream,outStream,order);
119 
120 } // checkGradient
121 
122 
123 
124 template <class Real>
125 std::vector<std::vector<Real> > Objective<Real>::checkGradient( const Vector<Real> &x,
126  const Vector<Real> &g,
127  const Vector<Real> &d,
128  const std::vector<Real> &steps,
129  const bool printToStream,
130  std::ostream & outStream,
131  const int order ) {
132 
133  TEUCHOS_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
134  "Error: finite difference order must be 1,2,3, or 4" );
135 
138 
139  Real tol = std::sqrt(ROL_EPSILON);
140 
141  int numSteps = steps.size();
142  int numVals = 4;
143  std::vector<Real> tmp(numVals);
144  std::vector<std::vector<Real> > gCheck(numSteps, tmp);
145 
146  // Save the format state of the original outStream.
147  Teuchos::oblackholestream oldFormatState;
148  oldFormatState.copyfmt(outStream);
149 
150  // Evaluate objective value at x.
151  this->update(x);
152  Real val = this->value(x,tol);
153 
154  // Compute gradient at x.
155  Teuchos::RCP<Vector<Real> > gtmp = g.clone();
156  this->gradient(*gtmp, x, tol);
157  Real dtg = d.dot(gtmp->dual());
158 
159  // Temporary vectors.
160  Teuchos::RCP<Vector<Real> > xnew = x.clone();
161 
162  for (int i=0; i<numSteps; i++) {
163 
164  Real eta = steps[i];
165 
166  xnew->set(x);
167 
168  // Compute gradient, finite-difference gradient, and absolute error.
169  gCheck[i][0] = eta;
170  gCheck[i][1] = dtg;
171 
172  gCheck[i][2] = weights[order-1][0] * val;
173 
174  for(int j=0; j<order; ++j) {
175  // Evaluate at x <- x+eta*c_i*d.
176  xnew->axpy(eta*shifts[order-1][j], d);
177 
178  // Only evaluate at shifts where the weight is nonzero
179  if( weights[order-1][j+1] != 0 ) {
180  this->update(*xnew);
181  gCheck[i][2] += weights[order-1][j+1] * this->value(*xnew,tol);
182  }
183  }
184 
185  gCheck[i][2] /= eta;
186 
187  gCheck[i][3] = std::abs(gCheck[i][2] - gCheck[i][1]);
188 
189  if (printToStream) {
190  if (i==0) {
191  outStream << std::right
192  << std::setw(20) << "Step size"
193  << std::setw(20) << "grad'*dir"
194  << std::setw(20) << "FD approx"
195  << std::setw(20) << "abs error"
196  << "\n";
197  }
198  outStream << std::scientific << std::setprecision(11) << std::right
199  << std::setw(20) << gCheck[i][0]
200  << std::setw(20) << gCheck[i][1]
201  << std::setw(20) << gCheck[i][2]
202  << std::setw(20) << gCheck[i][3]
203  << "\n";
204  }
205 
206  }
207 
208  // Reset format state of outStream.
209  outStream.copyfmt(oldFormatState);
210 
211  return gCheck;
212 } // checkGradient
213 
214 
215 
216 
217 
218 
219 
220 
221 template <class Real>
222 std::vector<std::vector<Real> > Objective<Real>::checkHessVec( const Vector<Real> &x,
223  const Vector<Real> &hv,
224  const Vector<Real> &v,
225  const bool printToStream,
226  std::ostream & outStream,
227  const int numSteps,
228  const int order ) {
229  std::vector<Real> steps(numSteps);
230  for(int i=0;i<numSteps;++i) {
231  steps[i] = pow(10,-i);
232  }
233 
234  return checkHessVec(x,hv,v,steps,printToStream,outStream,order);
235 } // checkHessVec
236 
237 
238 
239 template <class Real>
240 std::vector<std::vector<Real> > Objective<Real>::checkHessVec( const Vector<Real> &x,
241  const Vector<Real> &hv,
242  const Vector<Real> &v,
243  const std::vector<Real> &steps,
244  const bool printToStream,
245  std::ostream & outStream,
246  const int order ) {
247 
248  TEUCHOS_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
249  "Error: finite difference order must be 1,2,3, or 4" );
250 
253 
254 
255  Real tol = std::sqrt(ROL_EPSILON);
256 
257  int numSteps = steps.size();
258  int numVals = 4;
259  std::vector<Real> tmp(numVals);
260  std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
261 
262  // Save the format state of the original outStream.
263  Teuchos::oblackholestream oldFormatState;
264  oldFormatState.copyfmt(outStream);
265 
266  // Compute gradient at x.
267  Teuchos::RCP<Vector<Real> > g = hv.clone();
268  this->update(x);
269  this->gradient(*g, x, tol);
270 
271  // Compute (Hessian at x) times (vector v).
272  Teuchos::RCP<Vector<Real> > Hv = hv.clone();
273  this->hessVec(*Hv, v, x, tol);
274  Real normHv = Hv->norm();
275 
276  // Temporary vectors.
277  Teuchos::RCP<Vector<Real> > gdif = hv.clone();
278  Teuchos::RCP<Vector<Real> > gnew = hv.clone();
279  Teuchos::RCP<Vector<Real> > xnew = x.clone();
280 
281  for (int i=0; i<numSteps; i++) {
282 
283  Real eta = steps[i];
284 
285  // Evaluate objective value at x+eta*d.
286  xnew->set(x);
287 
288  gdif->set(*g);
289  gdif->scale(weights[order-1][0]);
290 
291  for(int j=0; j<order; ++j) {
292 
293  // Evaluate at x <- x+eta*c_i*d.
294  xnew->axpy(eta*shifts[order-1][j], v);
295 
296  // Only evaluate at shifts where the weight is nonzero
297  if( weights[order-1][j+1] != 0 ) {
298  this->update(*xnew);
299  this->gradient(*gnew, *xnew, tol);
300  gdif->axpy(weights[order-1][j+1],*gnew);
301  }
302 
303  }
304 
305  gdif->scale(1.0/eta);
306 
307  // Compute norms of hessvec, finite-difference hessvec, and error.
308  hvCheck[i][0] = eta;
309  hvCheck[i][1] = normHv;
310  hvCheck[i][2] = gdif->norm();
311  gdif->axpy(-1.0, *Hv);
312  hvCheck[i][3] = gdif->norm();
313 
314  if (printToStream) {
315  if (i==0) {
316  outStream << std::right
317  << std::setw(20) << "Step size"
318  << std::setw(20) << "norm(Hess*vec)"
319  << std::setw(20) << "norm(FD approx)"
320  << std::setw(20) << "norm(abs error)"
321  << "\n";
322  }
323  outStream << std::scientific << std::setprecision(11) << std::right
324  << std::setw(20) << hvCheck[i][0]
325  << std::setw(20) << hvCheck[i][1]
326  << std::setw(20) << hvCheck[i][2]
327  << std::setw(20) << hvCheck[i][3]
328  << "\n";
329  }
330 
331  }
332 
333  // Reset format state of outStream.
334  outStream.copyfmt(oldFormatState);
335 
336  return hvCheck;
337 } // checkHessVec
338 
339 
340 
341 template<class Real>
342 std::vector<Real> Objective<Real>::checkHessSym( const Vector<Real> &x,
343  const Vector<Real> &hv,
344  const Vector<Real> &v,
345  const Vector<Real> &w,
346  const bool printToStream,
347  std::ostream & outStream ) {
348 
349  Real tol = std::sqrt(ROL_EPSILON);
350 
351  // Compute (Hessian at x) times (vector v).
352  Teuchos::RCP<Vector<Real> > h = hv.clone();
353  this->hessVec(*h, v, x, tol);
354  Real wHv = w.dot(h->dual());
355 
356  this->hessVec(*h, w, x, tol);
357  Real vHw = v.dot(h->dual());
358 
359  std::vector<Real> hsymCheck(3, 0);
360 
361  hsymCheck[0] = wHv;
362  hsymCheck[1] = vHw;
363  hsymCheck[2] = std::abs(vHw-wHv);
364 
365  // Save the format state of the original outStream.
366  Teuchos::oblackholestream oldFormatState;
367  oldFormatState.copyfmt(outStream);
368 
369  if (printToStream) {
370  outStream << std::right
371  << std::setw(20) << "<w, H(x)v>"
372  << std::setw(20) << "<v, H(x)w>"
373  << std::setw(20) << "abs error"
374  << "\n";
375  outStream << std::scientific << std::setprecision(11) << std::right
376  << std::setw(20) << hsymCheck[0]
377  << std::setw(20) << hsymCheck[1]
378  << std::setw(20) << hsymCheck[2]
379  << "\n";
380  }
381 
382  // Reset format state of outStream.
383  outStream.copyfmt(oldFormatState);
384 
385  return hsymCheck;
386 
387 } // checkHessSym
388 
389 
390 
391 } // namespace ROL
392 
393 #endif
virtual void scale(const Real alpha)=0
Compute where .
virtual int dimension() const
Return dimension of the vector space.
Definition: ROL_Vector.hpp:181
const double weights[4][5]
Definition: ROL_Types.hpp:1001
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:141
virtual Real dirDeriv(const Vector< Real > &x, const Vector< Real > &d, Real &tol)
Compute directional derivative.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:155
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
virtual Real dot(const Vector &x) const =0
Compute where .
virtual std::vector< std::vector< Real > > checkGradient(const Vector< Real > &x, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference gradient check.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
virtual std::vector< std::vector< Real > > checkHessVec(const Vector< Real > &x, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference Hessian-applied-to-vector check.
virtual Teuchos::RCP< Vector > basis(const int i) const
Return i-th basis vector.
Definition: ROL_Vector.hpp:170
virtual Real norm() const =0
Returns where .
virtual std::vector< Real > checkHessSym(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout)
Hessian symmetry check.
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:115