ROL
ROL_SimulatedObjective.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_SIMULATED_OBJECTIVE_H
45 #define ROL_SIMULATED_OBJECTIVE_H
46 
47 #include "ROL_SimulatedVector.hpp"
48 #include "ROL_Objective_SimOpt.hpp"
49 
50 namespace ROL {
51 
52 template <class Real>
53 class SimulatedObjective : public Objective<Real> {
54 private:
55  const ROL::Ptr<SampleGenerator<Real> > sampler_;
56  const ROL::Ptr<Objective_SimOpt<Real> > pobj_;
57 
58 public:
59 
60  virtual ~SimulatedObjective() {}
61 
62  SimulatedObjective(const ROL::Ptr<SampleGenerator<Real> > & sampler,
63  const ROL::Ptr<Objective_SimOpt<Real> > & pobj)
64  : sampler_(sampler), pobj_(pobj) {}
65 
66  void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
67  pobj_->update(x,flag,iter);
68  }
69  void update( const Vector<Real> &x, UpdateType type, int iter = -1 ) {
70  pobj_->update(x,type,iter);
71  }
72 
73 
74  Real value(const Vector<Real> &x,
75  Real &tol) {
76  const Vector_SimOpt<Real> &uz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
77  ROL::Ptr<const Vector<Real> > uptr = uz.get_1();
78  ROL::Ptr<const Vector<Real> > zptr = uz.get_2();
79  const SimulatedVector<Real> &pu = dynamic_cast<const SimulatedVector<Real>&>(*uptr);
80 
81  std::vector<Real> param;
82  Real weight(0);
83  Real val = 0;
84  Real tmpval = 0;
85  Real tmpsum = 0;
86  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pu.numVectors(); ++i) {
87  param = sampler_->getMyPoint(static_cast<int>(i));
88  weight = sampler_->getMyWeight(static_cast<int>(i));
89  pobj_->setParameter(param);
90  pobj_->update(*(pu.get(i)), *zptr);
91  tmpval = pobj_->value(*(pu.get(i)), *zptr, tol);
92  tmpsum += tmpval*weight;
93  }
94  sampler_->sumAll(&tmpsum, &val, 1);
95  return val;
96  }
97 
98  virtual void gradient(Vector<Real> &g,
99  const Vector<Real> &x,
100  Real &tol) {
101  g.zero();
102  // split x
103  const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
104  ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
105  ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
106  const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
107  // split g
108  Vector_SimOpt<Real> &guz = dynamic_cast<Vector_SimOpt<Real>&>(g);
109  ROL::Ptr<Vector<Real> > guptr = guz.get_1();
110  ROL::Ptr<Vector<Real> > gzptr = guz.get_2();
111  SimulatedVector<Real> &pgu = dynamic_cast<SimulatedVector<Real>&>(*guptr);
112 
113  std::vector<Real> param;
114  Real weight(0);
115  ROL::Ptr<Vector<Real> > tmp1 = gzptr->clone();
116  ROL::Ptr<Vector<Real> > tmp2 = gzptr->clone();
117  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pgu.numVectors(); ++i) {
118  param = sampler_->getMyPoint(static_cast<int>(i));
119  weight = sampler_->getMyWeight(static_cast<int>(i));
120  pobj_->setParameter(param);
121  Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
122  Vector_SimOpt<Real> gi(pgu.get(i), tmp1);
123  pobj_->update(xi);
124  pobj_->gradient(gi, xi, tol);
125  gi.scale(weight);
126  tmp2->plus(*tmp1);
127  }
128  sampler_->sumAll(*tmp2, *gzptr);
129 
130  }
131 
132 
133  virtual void hessVec(Vector<Real> &hv,
134  const Vector<Real> &v,
135  const Vector<Real> &x,
136  Real &tol) {
137  hv.zero();
138  // split x
139  const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
140  ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
141  ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
142  const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
143  // split v
144  const Vector_SimOpt<Real> &vuz = dynamic_cast<const Vector_SimOpt<Real>&>(v);
145  ROL::Ptr<const Vector<Real> > vuptr = vuz.get_1();
146  ROL::Ptr<const Vector<Real> > vzptr = vuz.get_2();
147  const SimulatedVector<Real> &pvu = dynamic_cast<const SimulatedVector<Real>&>(*vuptr);
148  // split hv
149  Vector_SimOpt<Real> &hvuz = dynamic_cast<Vector_SimOpt<Real>&>(hv);
150  ROL::Ptr<Vector<Real> > hvuptr = hvuz.get_1();
151  ROL::Ptr<Vector<Real> > hvzptr = hvuz.get_2();
152  SimulatedVector<Real> &phvu = dynamic_cast<SimulatedVector<Real>&>(*hvuptr);
153 
154  std::vector<Real> param;
155  Real weight(0);
156  ROL::Ptr<Vector<Real> > tmp1 = hvzptr->clone();
157  ROL::Ptr<Vector<Real> > tmp2 = hvzptr->clone();
158  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<phvu.numVectors(); ++i) {
159  param = sampler_->getMyPoint(static_cast<int>(i));
160  weight = sampler_->getMyWeight(static_cast<int>(i));
161  pobj_->setParameter(param);
162  Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
163  Vector_SimOpt<Real> vi(ROL::constPtrCast<Vector<Real> >(pvu.get(i)), ROL::constPtrCast<Vector<Real> >(vzptr));
164  Vector_SimOpt<Real> hvi(phvu.get(i), tmp1);
165  pobj_->update(xi);
166  pobj_->hessVec(hvi, vi, xi, tol);
167  hvi.scale(weight);
168  tmp2->plus(*tmp1);
169  }
170  sampler_->sumAll(*tmp2, *hvzptr);
171  }
172 
173 }; // class SimulatedObjective
174 
175 } // namespace ROL
176 
177 #endif
Provides the interface to evaluate objective functions.
Provides the interface to evaluate simulation-based objective functions.
typename PV< Real >::size_type size_type
const ROL::Ptr< SampleGenerator< Real > > sampler_
ROL::Ptr< const Vector< Real > > get_2() const
Defines the linear algebra or vector space interface for simulation-based optimization.
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
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 the linear algebra of a vector space on a generic partitioned vector where the individual vec...
ROL::Ptr< const Vector< Real > > get(size_type i) const
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
size_type numVectors() const
const ROL::Ptr< Objective_SimOpt< Real > > pobj_
SimulatedObjective(const ROL::Ptr< SampleGenerator< Real > > &sampler, const ROL::Ptr< Objective_SimOpt< Real > > &pobj)
void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
ROL::Ptr< const Vector< Real > > get_1() const