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 
68  Real value(const Vector<Real> &x,
69  Real &tol) {
70  const Vector_SimOpt<Real> &uz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
71  ROL::Ptr<const Vector<Real> > uptr = uz.get_1();
72  ROL::Ptr<const Vector<Real> > zptr = uz.get_2();
73  const SimulatedVector<Real> &pu = dynamic_cast<const SimulatedVector<Real>&>(*uptr);
74 
75  std::vector<Real> param;
76  Real weight(0);
77  Real val = 0;
78  Real tmpval = 0;
79  Real tmpsum = 0;
80  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pu.numVectors(); ++i) {
81  param = sampler_->getMyPoint(static_cast<int>(i));
82  weight = sampler_->getMyWeight(static_cast<int>(i));
83  pobj_->setParameter(param);
84  pobj_->update(*(pu.get(i)), *zptr);
85  tmpval = pobj_->value(*(pu.get(i)), *zptr, tol);
86  tmpsum += tmpval*weight;
87  }
88  sampler_->sumAll(&tmpsum, &val, 1);
89  return val;
90  }
91 
92  virtual void gradient(Vector<Real> &g,
93  const Vector<Real> &x,
94  Real &tol) {
95  g.zero();
96  // split x
97  const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
98  ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
99  ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
100  const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
101  // split g
102  Vector_SimOpt<Real> &guz = dynamic_cast<Vector_SimOpt<Real>&>(g);
103  ROL::Ptr<Vector<Real> > guptr = guz.get_1();
104  ROL::Ptr<Vector<Real> > gzptr = guz.get_2();
105  SimulatedVector<Real> &pgu = dynamic_cast<SimulatedVector<Real>&>(*guptr);
106 
107  std::vector<Real> param;
108  Real weight(0);
109  ROL::Ptr<Vector<Real> > tmp1 = gzptr->clone();
110  ROL::Ptr<Vector<Real> > tmp2 = gzptr->clone();
111  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pgu.numVectors(); ++i) {
112  param = sampler_->getMyPoint(static_cast<int>(i));
113  weight = sampler_->getMyWeight(static_cast<int>(i));
114  pobj_->setParameter(param);
115  Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
116  Vector_SimOpt<Real> gi(pgu.get(i), tmp1);
117  pobj_->update(xi);
118  pobj_->gradient(gi, xi, tol);
119  gi.scale(weight);
120  tmp2->plus(*tmp1);
121  }
122  sampler_->sumAll(*tmp2, *gzptr);
123 
124  }
125 
126 
127  virtual void hessVec(Vector<Real> &hv,
128  const Vector<Real> &v,
129  const Vector<Real> &x,
130  Real &tol) {
131  hv.zero();
132  // split x
133  const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
134  ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
135  ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
136  const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
137  // split v
138  const Vector_SimOpt<Real> &vuz = dynamic_cast<const Vector_SimOpt<Real>&>(v);
139  ROL::Ptr<const Vector<Real> > vuptr = vuz.get_1();
140  ROL::Ptr<const Vector<Real> > vzptr = vuz.get_2();
141  const SimulatedVector<Real> &pvu = dynamic_cast<const SimulatedVector<Real>&>(*vuptr);
142  // split hv
143  Vector_SimOpt<Real> &hvuz = dynamic_cast<Vector_SimOpt<Real>&>(hv);
144  ROL::Ptr<Vector<Real> > hvuptr = hvuz.get_1();
145  ROL::Ptr<Vector<Real> > hvzptr = hvuz.get_2();
146  SimulatedVector<Real> &phvu = dynamic_cast<SimulatedVector<Real>&>(*hvuptr);
147 
148  std::vector<Real> param;
149  Real weight(0);
150  ROL::Ptr<Vector<Real> > tmp1 = hvzptr->clone();
151  ROL::Ptr<Vector<Real> > tmp2 = hvzptr->clone();
152  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<phvu.numVectors(); ++i) {
153  param = sampler_->getMyPoint(static_cast<int>(i));
154  weight = sampler_->getMyWeight(static_cast<int>(i));
155  pobj_->setParameter(param);
156  Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
157  Vector_SimOpt<Real> vi(ROL::constPtrCast<Vector<Real> >(pvu.get(i)), ROL::constPtrCast<Vector<Real> >(vzptr));
158  Vector_SimOpt<Real> hvi(phvu.get(i), tmp1);
159  pobj_->update(xi);
160  pobj_->hessVec(hvi, vi, xi, tol);
161  hvi.scale(weight);
162  tmp2->plus(*tmp1);
163  }
164  sampler_->sumAll(*tmp2, *hvzptr);
165  }
166 
167 }; // class SimulatedObjective
168 
169 } // namespace ROL
170 
171 #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)
Real value(const Vector< Real > &x, Real &tol)
Compute value.
ROL::Ptr< const Vector< Real > > get_1() const