ROL
example_01b.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 
54 #include "ROL_StdVector.hpp"
55 #include "Sacado.hpp"
56 
57 using namespace ROL;
58 
59 template<class ScalarT>
61 
62  public:
63 
64  ScalarT eval(const std::vector<ScalarT> &x);
65 
66 };
67 
78 template<class ScalarT>
79 ScalarT FunctionZakharov<ScalarT>::eval(const std::vector<ScalarT> & x) {
80 
81  ScalarT xdotx = 0;
82  ScalarT kdotx = 0;
83  ScalarT J = 0;
84 
85  // Compute dot products
86  for(unsigned i=0; i<x.size(); ++i) {
87  xdotx += pow(x[i],2); // (k,x)
88  kdotx += double(i+1)*x[i]; // (x,x)
89  }
90 
91  // Sum terms in objective function
92  J = xdotx + pow(kdotx,2)/4.0 + pow(kdotx,4)/16.0;
93 
94  return J;
95 }
96 
97 
98 template<class Real>
99 class Zakharov_Sacado_Objective : public Objective<Real> {
100 
101  typedef Sacado::Fad::DFad<Real> GradType;
102  typedef Sacado::Fad::SFad<Real,1> DirDerivType;
103  typedef Sacado::Fad::DFad<DirDerivType> HessVecType;
104  // In C++11, we could use template typedefs:
105  // template <typename T> using GradTypeT = Sacado::Fad::DFad<T>;
106  // typedef Sacado::Fad::SFad<Real,1> DirDerivType;
107  // typedef GradTypeT<Real> GradType;
108  // typedef GradTypeT<DirDerivType> HessVecType;
109 
110  private:
111 
115 
116  public:
117 
119 
120  /* Evaluate the objective function at x */
121  Real value( const Vector<Real> &x, Real &tol ) {
122  Teuchos::RCP<const std::vector<Real> > xp =
123  (Teuchos::dyn_cast<const StdVector<Real> >(x)).getVector();
124  return zfunc_.eval(*xp);
125  }
126 
127  /* Evaluate the gradient at x */
128  void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
129  Teuchos::RCP<std::vector<Real> > gp =
130  (Teuchos::dyn_cast<StdVector<Real> >(g)).getVector();
131  Teuchos::RCP<const std::vector<Real> > xp =
132  (Teuchos::dyn_cast<const StdVector<Real> >(x)).getVector();
133 
134  unsigned n = xp->size();
135 
136  std::vector<GradType> x_grad(n);
137 
138  for(unsigned i=0; i<n; ++i) {
139  x_grad[i] = (*xp)[i]; // Set values x(i).
140  x_grad[i].diff(i,n); // Choose canonical directions.
141  }
142 
143  GradType J_grad = zfuncGrad_.eval(x_grad);
144 
145  for(unsigned i=0; i<n; ++i) {
146  (*gp)[i] = J_grad.dx(i);
147  }
148 
149  }
150 
151  /* Compute the action of the Hessian evaluated at x on a vector v */
152  void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
153  Teuchos::RCP<std::vector<Real> > hvp =
154  (Teuchos::dyn_cast<StdVector<Real> >(hv)).getVector();
155  Teuchos::RCP<const std::vector<Real> > vp =
156  (Teuchos::dyn_cast<const StdVector<Real> >(v)).getVector();
157  Teuchos::RCP<const std::vector<Real> > xp =
158  (Teuchos::dyn_cast<const StdVector<Real> >(x)).getVector();
159 
160  unsigned n = xp->size();
161 
162  std::vector<HessVecType> x_hessvec(n);
163 
164  for(unsigned i=0; i<n; ++i) {
165  DirDerivType tmp(1,(*xp)[i]); // Set values x(i).
166  tmp.fastAccessDx(0)= (*vp)[i]; // Set direction values v(i).
167  x_hessvec[i] = tmp; // Use tmp to define hessvec-able x.
168  x_hessvec[i].diff(i,n); // Choose directions.
169  }
170 
171  // Compute Hessian-vector product (and other currently irrelevant things).
172  HessVecType J_hessvec = zfuncHessVec_.eval(x_hessvec);
173 
174  for(unsigned i=0; i<n; ++i) {
175  (*hvp)[i] = (J_hessvec.dx(i)).fastAccessDx(0);
176  // hessvec = get gradient get dir deriv
177  }
178  }
179 };
180 
Provides the interface to evaluate objective functions.
FunctionZakharov< HessVecType > zfuncHessVec_
ScalarT eval(const std::vector< ScalarT > &x)
A Sacado-accessible version of the Zakharov function to differentiate Where .
Definition: example_01b.hpp:79
FunctionZakharov< Real > zfunc_
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
FunctionZakharov< GradType > zfuncGrad_
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
Provides the std::vector implementation of the ROL::Vector interface.
Sacado::Fad::SFad< Real, 1 > DirDerivType
Sacado::Fad::DFad< Real > GradType
Real value(const Vector< Real > &x, Real &tol)
Compute value.
Sacado::Fad::DFad< DirDerivType > HessVecType
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.