ROL
ROL_SimulatedObjective.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_SIMULATED_OBJECTIVE_H
11#define ROL_SIMULATED_OBJECTIVE_H
12
15
16namespace ROL {
17
18template <class Real>
19class SimulatedObjective : public Objective<Real> {
20private:
21 const ROL::Ptr<SampleGenerator<Real> > sampler_;
22 const ROL::Ptr<Objective_SimOpt<Real> > pobj_;
23
24public:
25
27
28 SimulatedObjective(const ROL::Ptr<SampleGenerator<Real> > & sampler,
29 const ROL::Ptr<Objective_SimOpt<Real> > & pobj)
30 : sampler_(sampler), pobj_(pobj) {}
31
32 void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
33 pobj_->update(x,flag,iter);
34 }
35 void update( const Vector<Real> &x, UpdateType type, int iter = -1 ) {
36 pobj_->update(x,type,iter);
37 }
38
39
40 Real value(const Vector<Real> &x,
41 Real &tol) {
42 const Vector_SimOpt<Real> &uz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
43 ROL::Ptr<const Vector<Real> > uptr = uz.get_1();
44 ROL::Ptr<const Vector<Real> > zptr = uz.get_2();
45 const SimulatedVector<Real> &pu = dynamic_cast<const SimulatedVector<Real>&>(*uptr);
46
47 std::vector<Real> param;
48 Real weight(0);
49 Real val = 0;
50 Real tmpval = 0;
51 Real tmpsum = 0;
52 for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pu.numVectors(); ++i) {
53 param = sampler_->getMyPoint(static_cast<int>(i));
54 weight = sampler_->getMyWeight(static_cast<int>(i));
55 pobj_->setParameter(param);
56 pobj_->update(*(pu.get(i)), *zptr);
57 tmpval = pobj_->value(*(pu.get(i)), *zptr, tol);
58 tmpsum += tmpval*weight;
59 }
60 sampler_->sumAll(&tmpsum, &val, 1);
61 return val;
62 }
63
64 virtual void gradient(Vector<Real> &g,
65 const Vector<Real> &x,
66 Real &tol) {
67 g.zero();
68 // split x
69 const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
70 ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
71 ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
72 const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
73 // split g
74 Vector_SimOpt<Real> &guz = dynamic_cast<Vector_SimOpt<Real>&>(g);
75 ROL::Ptr<Vector<Real> > guptr = guz.get_1();
76 ROL::Ptr<Vector<Real> > gzptr = guz.get_2();
77 SimulatedVector<Real> &pgu = dynamic_cast<SimulatedVector<Real>&>(*guptr);
78
79 std::vector<Real> param;
80 Real weight(0);
81 ROL::Ptr<Vector<Real> > tmp1 = gzptr->clone();
82 ROL::Ptr<Vector<Real> > tmp2 = gzptr->clone();
83 for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pgu.numVectors(); ++i) {
84 param = sampler_->getMyPoint(static_cast<int>(i));
85 weight = sampler_->getMyWeight(static_cast<int>(i));
86 pobj_->setParameter(param);
87 Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
88 Vector_SimOpt<Real> gi(pgu.get(i), tmp1);
89 pobj_->update(xi);
90 pobj_->gradient(gi, xi, tol);
91 gi.scale(weight);
92 tmp2->plus(*tmp1);
93 }
94 sampler_->sumAll(*tmp2, *gzptr);
95
96 }
97
98
99 virtual void hessVec(Vector<Real> &hv,
100 const Vector<Real> &v,
101 const Vector<Real> &x,
102 Real &tol) {
103 hv.zero();
104 // split x
105 const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
106 ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
107 ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
108 const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
109 // split v
110 const Vector_SimOpt<Real> &vuz = dynamic_cast<const Vector_SimOpt<Real>&>(v);
111 ROL::Ptr<const Vector<Real> > vuptr = vuz.get_1();
112 ROL::Ptr<const Vector<Real> > vzptr = vuz.get_2();
113 const SimulatedVector<Real> &pvu = dynamic_cast<const SimulatedVector<Real>&>(*vuptr);
114 // split hv
115 Vector_SimOpt<Real> &hvuz = dynamic_cast<Vector_SimOpt<Real>&>(hv);
116 ROL::Ptr<Vector<Real> > hvuptr = hvuz.get_1();
117 ROL::Ptr<Vector<Real> > hvzptr = hvuz.get_2();
118 SimulatedVector<Real> &phvu = dynamic_cast<SimulatedVector<Real>&>(*hvuptr);
119
120 std::vector<Real> param;
121 Real weight(0);
122 ROL::Ptr<Vector<Real> > tmp1 = hvzptr->clone();
123 ROL::Ptr<Vector<Real> > tmp2 = hvzptr->clone();
124 for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<phvu.numVectors(); ++i) {
125 param = sampler_->getMyPoint(static_cast<int>(i));
126 weight = sampler_->getMyWeight(static_cast<int>(i));
127 pobj_->setParameter(param);
128 Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
129 Vector_SimOpt<Real> vi(ROL::constPtrCast<Vector<Real> >(pvu.get(i)), ROL::constPtrCast<Vector<Real> >(vzptr));
130 Vector_SimOpt<Real> hvi(phvu.get(i), tmp1);
131 pobj_->update(xi);
132 pobj_->hessVec(hvi, vi, xi, tol);
133 hvi.scale(weight);
134 tmp2->plus(*tmp1);
135 }
136 sampler_->sumAll(*tmp2, *hvzptr);
137 }
138
139}; // class SimulatedObjective
140
141} // namespace ROL
142
143#endif
typename PV< Real >::size_type size_type
Provides the interface to evaluate simulation-based objective functions.
Provides the interface to evaluate objective functions.
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.
const ROL::Ptr< SampleGenerator< Real > > sampler_
const ROL::Ptr< Objective_SimOpt< Real > > pobj_
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
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
Defines the linear algebra or vector space interface for simulation-based optimization.
ROL::Ptr< const Vector< Real > > get_2() const
ROL::Ptr< const Vector< Real > > get_1() const
void scale(const Real alpha)
Compute where .
Defines the linear algebra or vector space interface.
virtual void zero()
Set to zero vector.