ROL
function/test_08.cpp
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 #include "ROL_Stream.hpp"
45 #include "Teuchos_GlobalMPISession.hpp"
46 
47 #include "ROL_StdVector.hpp"
48 #include "ROL_StdObjective.hpp"
49 #include "ROL_Objective_SimOpt.hpp"
51 
52 typedef double RealT;
53 
54 template<class Real>
56 public:
57  Real value( const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
58  ROL::Ptr<const std::vector<Real> > up
59  = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
60  ROL::Ptr<const std::vector<Real> > zp
61  = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
62  Real half(0.5), quadu(0), quadz(0);
63  unsigned usize = up->size();
64  for ( unsigned i = 0; i < usize; i++ ) {
65  quadu += (*up)[i]*(*up)[i];
66  }
67  unsigned zsize = zp->size();
68  for ( unsigned i = 0; i < zsize; i++ ) {
69  quadz += (*zp)[i]*(*zp)[i];
70  }
71  return half*(quadu + quadz);
72  }
73 
74  void gradient_1( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
75  ROL::Ptr<std::vector<Real> > gp
76  = dynamic_cast<ROL::StdVector<Real>&>(g).getVector();
77  ROL::Ptr<const std::vector<Real> > up
78  = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
79  gp->assign(up->begin(),up->end());
80  }
81 
82  void gradient_2( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
83  ROL::Ptr<std::vector<Real> > gp
84  = dynamic_cast<ROL::StdVector<Real>&>(g).getVector();
85  ROL::Ptr<const std::vector<Real> > zp
86  = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
87  gp->assign(zp->begin(),zp->end());
88  }
89 
91  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
92  ROL::Ptr<std::vector<Real> > hvp
93  = dynamic_cast<ROL::StdVector<Real>&>(hv).getVector();
94  ROL::Ptr<const std::vector<Real> > vp
95  = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
96  hvp->assign(vp->begin(),vp->end());
97  }
98 
100  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
101  hv.zero();
102  }
103 
105  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
106  hv.zero();
107  }
108 
110  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
111  ROL::Ptr<std::vector<Real> > hvp
112  = dynamic_cast<ROL::StdVector<Real>&>(hv).getVector();
113  ROL::Ptr<const std::vector<Real> > vp
114  = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
115  hvp->assign(vp->begin(),vp->end());
116  }
117 };
118 
119 template<class Real>
121 public:
122  Real value( const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
123  ROL::Ptr<const std::vector<Real> > up
124  = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
125  ROL::Ptr<const std::vector<Real> > zp
126  = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
127  Real linu(0), linz(0);
128  unsigned usize = up->size();
129  for ( unsigned i = 0; i < usize; i++ ) {
130  linu += (*up)[i];
131  }
132  unsigned zsize = zp->size();
133  for ( unsigned i = 0; i < zsize; i++ ) {
134  linz += (*zp)[i];
135  }
136  return linu + linz;
137  }
138 
139  void gradient_1( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
140  ROL::Ptr<std::vector<Real> > gp
141  = dynamic_cast<ROL::StdVector<Real>&>(g).getVector();
142  ROL::Ptr<const std::vector<Real> > up
143  = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
144  gp->assign(up->size(),1);
145  }
146 
147  void gradient_2( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
148  ROL::Ptr<std::vector<Real> > gp
149  = dynamic_cast<ROL::StdVector<Real>&>(g).getVector();
150  ROL::Ptr<const std::vector<Real> > zp
151  = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
152  gp->assign(zp->size(),1);
153  }
154 
156  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
157  hv.zero();
158  }
159 
161  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
162  hv.zero();
163  }
164 
166  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
167  hv.zero();
168  }
169 
171  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
172  hv.zero();
173  }
174 };
175 
176 template<class Real>
178 public:
179  Real value( const std::vector<Real> &x, Real &tol ) {
180  return std::log(x[0]) * std::exp(x[1]);
181  }
182 
183  void gradient( std::vector<Real> &g, const std::vector<Real> &x, Real &tol ) {
184  g[0] = std::exp(x[1])/x[0];
185  g[1] = std::exp(x[1]) * std::log(x[0]);
186  }
187 
188  void hessVec( std::vector<Real> &hv, const std::vector<Real> &v, const std::vector<Real> &x, Real &tol ) {
189  Real H11 = -std::exp(x[1])/(x[0]*x[0]);
190  Real H12 = std::exp(x[1])/x[0];
191  Real H21 = std::exp(x[1])/x[0];
192  Real H22 = std::exp(x[1]) * std::log(x[0]);
193  hv[0] = H11*v[0] + H12*v[1];
194  hv[1] = H21*v[0] + H22*v[1];
195  }
196 };
197 
198 void setRandomVector(std::vector<RealT> &x) {
199  unsigned dim = x.size();
200  for ( unsigned i = 0; i < dim; i++ ) {
201  x[i] = (RealT)rand()/(RealT)RAND_MAX;
202  }
203 }
204 
205 int main(int argc, char* argv[]) {
206 
207  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
208 
209  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
210  int iprint = argc - 1;
211  ROL::Ptr<std::ostream> outStream;
212  ROL::nullstream bhs; // outputs nothing
213  if (iprint > 0)
214  outStream = ROL::makePtrFromRef(std::cout);
215  else
216  outStream = ROL::makePtrFromRef(bhs);
217 
218  int errorFlag = 0;
219 
220  try {
221  /**********************************************************************************************/
222  /************************* CONSTRUCT SOL COMPONENTS *******************************************/
223  /**********************************************************************************************/
224  // Build vectors
225  unsigned dim = 2;
226  ROL::Ptr<std::vector<RealT> > u_ptr = ROL::makePtr<std::vector<RealT>>(dim,0.0);
227  ROL::Ptr<ROL::Vector<RealT> > u = ROL::makePtr<ROL::StdVector<RealT>>(u_ptr);
228  setRandomVector(*u_ptr);
229  ROL::Ptr<std::vector<RealT> > z_ptr = ROL::makePtr<std::vector<RealT>>(dim,0.0);
230  ROL::Ptr<ROL::Vector<RealT> > z = ROL::makePtr<ROL::StdVector<RealT>>(z_ptr);
231  setRandomVector(*z_ptr);
232  ROL::Ptr<ROL::Vector<RealT> > x = ROL::makePtr<ROL::Vector_SimOpt<RealT>>(u,z);
233  ROL::Ptr<std::vector<RealT> > du_ptr = ROL::makePtr<std::vector<RealT>>(dim,0.0);
234  ROL::Ptr<ROL::Vector<RealT> > du = ROL::makePtr<ROL::StdVector<RealT>>(du_ptr);
235  setRandomVector(*du_ptr);
236  ROL::Ptr<std::vector<RealT> > dz_ptr = ROL::makePtr<std::vector<RealT>>(dim,0.0);
237  ROL::Ptr<ROL::Vector<RealT> > dz = ROL::makePtr<ROL::StdVector<RealT>>(dz_ptr);
238  setRandomVector(*dz_ptr);
239  ROL::Ptr<ROL::Vector<RealT> > d = ROL::makePtr<ROL::Vector_SimOpt<RealT>>(du,dz);
240  // Build objective function
241  std::vector<ROL::Ptr<ROL::Objective_SimOpt<RealT> > > vec_obj(2,ROL::nullPtr);
242  vec_obj[0] = ROL::makePtr<ObjectiveFunctionTest08_1<RealT>>();
243  vec_obj[1] = ROL::makePtr<ObjectiveFunctionTest08_2<RealT>>();
244  ROL::Ptr<ROL::StdObjective<RealT> > obj_scalarize
245  = ROL::makePtr<ObjectiveFunctionTest08_scalarize<RealT>>();
246  ROL::Ptr<ROL::CompositeObjective_SimOpt<RealT> > obj
247  = ROL::makePtr<ROL::CompositeObjective_SimOpt<RealT>>(vec_obj,obj_scalarize);
248  // Test parametrized objective functions
249  *outStream << "Check Derivatives of CompositeObjective_SimOpt\n";
250  obj->checkGradient(*x,*d,true,*outStream);
251  obj->checkHessVec(*x,*d,true,*outStream);
252  }
253  catch (std::logic_error& err) {
254  *outStream << err.what() << "\n";
255  errorFlag = -1000;
256  }; // end try
257 
258  if (errorFlag != 0)
259  std::cout << "End Result: TEST FAILED\n";
260  else
261  std::cout << "End Result: TEST PASSED\n";
262 
263  return 0;
264 }
Provides the interface to evaluate simulation-based objective functions.
void hessVec_12(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
void setRandomVector(std::vector< RealT > &x)
void gradient_2(ROL::Vector< Real > &g, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute gradient with respect to second component.
void hessVec_21(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
void hessVec(std::vector< Real > &hv, const std::vector< Real > &v, const std::vector< Real > &x, Real &tol)
void gradient_2(ROL::Vector< Real > &g, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute gradient with respect to second component.
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:167
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
Real value(const std::vector< Real > &x, Real &tol)
Specializes the ROL::Objective interface for objective functions that operate on ROL::StdVector&#39;s.
void gradient_1(ROL::Vector< Real > &g, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute gradient with respect to first component.
void hessVec_12(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
void hessVec_11(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
void hessVec_21(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
void hessVec_11(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
Real value(const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute value.
basic_nullstream< char, char_traits< char >> nullstream
Definition: ROL_Stream.hpp:72
int main(int argc, char *argv[])
void gradient_1(ROL::Vector< Real > &g, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute gradient with respect to first component.
void hessVec_22(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
void hessVec_22(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
constexpr auto dim
void gradient(std::vector< Real > &g, const std::vector< Real > &x, Real &tol)
Real value(const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute value.