ROL
ROL_MomentObjective.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_MOMENTOBJECTIVE_H
45#define ROL_MOMENTOBJECTIVE_H
46
47#include "ROL_Objective.hpp"
48#include "ROL_BatchManager.hpp"
49#include "ROL_Distribution.hpp"
50#include "ROL_SROMVector.hpp"
51#include "ROL_Types.hpp"
52#include <iostream>
53
54namespace ROL {
55
56template <class Real>
57class MomentObjective : public Objective<Real> {
58private:
59 std::vector<std::vector<std::pair<int, Real> > > moments_;
60 ROL::Ptr<BatchManager<Real> > bman_;
63 const bool optProb_;
64 const bool optAtom_;
65
66 Real momentValue(const int dim, const Real power, const Real moment,
67 const ProbabilityVector<Real> &prob,
68 const AtomVector<Real> &atom) const {
69 const int numSamples = prob.getNumMyAtoms();
70 Real val(0), xpt(0), xwt(0), sum(0), half(0.5), one(1), two(2);
71 for (int k = 0; k < numSamples; k++) {
72 xpt = (*atom.getAtom(k))[dim]; xwt = prob.getProbability(k);
73 val += xwt * ((power==one) ? xpt : std::pow(xpt,power));
74 }
75 bman_->sumAll(&val,&sum,1);
76 Real denom = ((std::abs(moment) < ROL_EPSILON<Real>()) ? one : moment);
77 return half*std::pow((sum-moment)/denom,two);
78 }
79
80 void momentGradient(std::vector<Real> &gradx, std::vector<Real> &gradp, Real &scale,
81 const int dim, const Real power, const Real moment,
82 const ProbabilityVector<Real> &prob,
83 const AtomVector<Real> &atom) const {
84 const int numSamples = prob.getNumMyAtoms();
85 gradx.resize(numSamples,0); gradp.resize(numSamples,0);
86 scale = 0;
87 Real xpt(0), xwt(0), xpow(0), psum(0), one(1), two(2);
88 for (int k = 0; k < numSamples; k++) {
89 xpt = (*atom.getAtom(k))[dim]; xwt = prob.getProbability(k);
90 xpow = ((power==one) ? one : ((power==two) ? xpt : std::pow(xpt,power-one)));
91 psum += xwt * xpow * xpt;
92 gradx[k] = xwt * xpow * power;
93 gradp[k] = xpow * xpt;
94 }
95 bman_->sumAll(&psum,&scale,1);
96 scale -= moment;
97 Real denom = ((std::abs(moment) < ROL_EPSILON<Real>()) ? one : moment);
98 scale /= std::pow(denom,two);
99 }
100
101 void momentHessVec(std::vector<Real> &hvx1, std::vector<Real> &hvx2, std::vector<Real> &hvx3,
102 std::vector<Real> &hvp1, std::vector<Real> &hvp2,
103 Real &scale1, Real &scale2, Real &scale3,
104 const int dim, const Real power, const Real moment,
105 const ProbabilityVector<Real> &prob,
106 const AtomVector<Real> &atom,
107 const ProbabilityVector<Real> &vprob,
108 const AtomVector<Real> &vatom) const {
109 const int numSamples = prob.getNumMyAtoms();
110 hvx1.resize(numSamples,0); hvx2.resize(numSamples,0); hvx3.resize(numSamples,0);
111 hvp1.resize(numSamples,0); hvp2.resize(numSamples,0);
112 scale1 = 0; scale2 = 0; scale3 = 0;
113 std::vector<Real> psum(3,0), scale(3,0);
114 Real xpt(0), xwt(0), vpt(0), vwt(0);
115 Real xpow0(0), xpow1(0), xpow2(0), zero(0), one(1), two(2), three(3);
116 for (int k = 0; k < numSamples; k++) {
117 xpt = (*atom.getAtom(k))[dim]; xwt = prob.getProbability(k);
118 vpt = (*vatom.getAtom(k))[dim]; vwt = vprob.getProbability(k);
119 xpow2 = ((power==one) ? zero : ((power==two) ? one : ((power==three) ? xpt :
120 std::pow(xpt,power-two))));
121 xpow1 = ((power==one) ? one : xpow2 * xpt);
122 xpow0 = xpow1 * xpt;
123 psum[0] += xwt * xpow1 * vpt;
124 psum[1] += xwt * xpow0;
125 psum[2] += vwt * xpow0;
126 hvx1[k] = power * xwt * xpow1;
127 hvx2[k] = power * (power-one) * xwt * xpow2 * vpt;
128 hvx3[k] = power * vwt * xpow1;
129 hvp1[k] = xpow0;
130 hvp2[k] = power * xpow1 * vpt;
131 }
132 bman_->sumAll(&psum[0],&scale[0],3);
133 Real denom = ((std::abs(moment) < ROL_EPSILON<Real>()) ? one : moment);
134 Real denom2 = denom*denom;
135 //const Real moment2 = std::pow(moment,2);
136 scale1 = scale[0] * power/denom2;
137 scale2 = (scale[1] - moment)/denom2 ;
138 scale3 = scale[2]/denom2;
139 }
140
141public:
142 MomentObjective(const std::vector<std::vector<std::pair<int, Real> > > &moments,
143 const ROL::Ptr<BatchManager<Real> > &bman,
144 const bool optProb = true, const bool optAtom = true)
145 : Objective<Real>(), moments_(moments), bman_(bman),
146 optProb_(optProb), optAtom_(optAtom) {
147 dimension_ = moments_.size();
148 numMoments_ = moments_[0].size();
149 }
150
151 MomentObjective(const std::vector<ROL::Ptr<Distribution<Real> > > &dist,
152 const std::vector<int> &order,
153 const ROL::Ptr<BatchManager<Real> > &bman,
154 const bool optProb = true, const bool optAtom = true)
155 : Objective<Real>(), bman_(bman), optProb_(optProb), optAtom_(optAtom) {
156 numMoments_ = order.size();
157 dimension_ = dist.size();
158 std::vector<std::pair<int,Real> > data(numMoments_);
159 moments_.clear(); moments_.resize(dimension_);
160 for (int d = 0; d < dimension_; d++) {
161 for (int i = 0; i < numMoments_; i++) {
162 data[i] = std::make_pair(order[i],dist[d]->moment(order[i]));
163 }
164 moments_[d].assign(data.begin(),data.end());
165 }
166 }
167
168 Real value( const Vector<Real> &x, Real &tol ) {
169 const SROMVector<Real> &ex = dynamic_cast<const SROMVector<Real>&>(x);
170 const ProbabilityVector<Real> &prob = *(ex.getProbabilityVector());
171 const AtomVector<Real> &atom = *(ex.getAtomVector());
172 Real val(0);
173 std::vector<std::pair<int, Real> > data;
174 for (int d = 0; d < dimension_; d++) {
175 data = moments_[d];
176 for (int m = 0; m < numMoments_; m++) {
177 val += momentValue(d,(Real)data[m].first,data[m].second,prob,atom);
178 }
179 }
180 return val;
181 }
182
183 void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
184 g.zero();
185 const SROMVector<Real> &ex = dynamic_cast<const SROMVector<Real>&>(x);
186 const ProbabilityVector<Real> &prob = *(ex.getProbabilityVector());
187 const AtomVector<Real> &atom = *(ex.getAtomVector());
188 int numSamples = prob.getNumMyAtoms();
189 std::vector<Real> gradx(numSamples,0), gradp(numSamples,0);
190 Real scale(0);
191 std::vector<std::pair<int, Real> > data;
192 std::vector<Real> val_wt(numSamples,0), tmp(dimension_,0);
193 std::vector<std::vector<Real> > val_pt(numSamples,tmp);
194 for (int d = 0; d < dimension_; d++) {
195 data = moments_[d];
196 for (int m = 0; m < numMoments_; m++) {
197 momentGradient(gradx,gradp,scale,d,(Real)data[m].first,data[m].second,prob,atom);
198 for (int k = 0; k < numSamples; k++) {
199 (val_pt[k])[d] += scale*gradx[k];
200 val_wt[k] += scale*gradp[k];
201 }
202 }
203 }
204 SROMVector<Real> &eg = dynamic_cast<SROMVector<Real>&>(g);
206 AtomVector<Real> &gatom = *(eg.getAtomVector());
207 for (int k = 0; k < numSamples; k++) {
208 if ( optProb_ ) {
209 gprob.setProbability(k,val_wt[k]);
210 }
211 if ( optAtom_ ) {
212 gatom.setAtom(k,val_pt[k]);
213 }
214 }
215 }
216
217 void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
218 hv.zero();
219 const SROMVector<Real> &ev = dynamic_cast<const SROMVector<Real>&>(v);
220 const ProbabilityVector<Real> &vprob = *(ev.getProbabilityVector());
221 const AtomVector<Real> &vatom = *(ev.getAtomVector());
222 const SROMVector<Real> &ex = dynamic_cast<const SROMVector<Real>&>(x);
223 const ProbabilityVector<Real> &prob = *(ex.getProbabilityVector());
224 const AtomVector<Real> &atom = *(ex.getAtomVector());
225 const int numSamples = prob.getNumMyAtoms();
226 std::vector<Real> hvx1(numSamples,0), hvx2(numSamples,0), hvx3(numSamples,0);
227 std::vector<Real> hvp1(numSamples,0), hvp2(numSamples,0);
228 Real scale1(0), scale2(0), scale3(0);
229 std::vector<std::pair<int, Real> > data;
230 std::vector<Real> val_wt(numSamples,0), tmp(dimension_,0);
231 std::vector<std::vector<Real> > val_pt(numSamples,tmp);
232 for (int d = 0; d < dimension_; d++) {
233 data = moments_[d];
234 for (int m = 0; m < numMoments_; m++) {
235 momentHessVec(hvx1,hvx2,hvx3,hvp1,hvp2,scale1,scale2,scale3,
236 d,(Real)data[m].first,data[m].second,prob,atom,vprob,vatom);
237 for (int k = 0; k < numSamples; k++) {
238 (val_pt[k])[d] += (scale1+scale3)*hvx1[k] + scale2*(hvx2[k]+hvx3[k]);
239 val_wt[k] += (scale1+scale3)*hvp1[k] + scale2*hvp2[k];
240 }
241 }
242 }
243 SROMVector<Real> &ehv = dynamic_cast<SROMVector<Real>&>(hv);
245 AtomVector<Real> &hatom = *(ehv.getAtomVector());
246 for (int k = 0; k < numSamples; k++) {
247 if ( optProb_ ) {
248 hprob.setProbability(k,val_wt[k]);
249 }
250 if ( optAtom_ ) {
251 hatom.setAtom(k,val_pt[k]);
252 }
253 }
254 }
255}; // class SROMObjective
256
257} // namespace ROL
258
259#endif
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Contains definitions of custom data types in ROL.
Provides the std::vector implementation of the ROL::Vector interface.
ROL::Ptr< const std::vector< Real > > getAtom(const int i) const
void setAtom(const int i, const std::vector< Real > &pt)
void momentHessVec(std::vector< Real > &hvx1, std::vector< Real > &hvx2, std::vector< Real > &hvx3, std::vector< Real > &hvp1, std::vector< Real > &hvp2, Real &scale1, Real &scale2, Real &scale3, const int dim, const Real power, const Real moment, const ProbabilityVector< Real > &prob, const AtomVector< Real > &atom, const ProbabilityVector< Real > &vprob, const AtomVector< Real > &vatom) const
MomentObjective(const std::vector< ROL::Ptr< Distribution< Real > > > &dist, const std::vector< int > &order, const ROL::Ptr< BatchManager< Real > > &bman, const bool optProb=true, const bool optAtom=true)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Real momentValue(const int dim, const Real power, const Real moment, const ProbabilityVector< Real > &prob, const AtomVector< Real > &atom) const
void momentGradient(std::vector< Real > &gradx, std::vector< Real > &gradp, Real &scale, const int dim, const Real power, const Real moment, const ProbabilityVector< Real > &prob, const AtomVector< Real > &atom) const
std::vector< std::vector< std::pair< int, Real > > > moments_
Real value(const Vector< Real > &x, Real &tol)
Compute value.
MomentObjective(const std::vector< std::vector< std::pair< int, Real > > > &moments, const ROL::Ptr< BatchManager< Real > > &bman, const bool optProb=true, const bool optAtom=true)
ROL::Ptr< BatchManager< Real > > bman_
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
Provides the interface to evaluate objective functions.
Provides the std::vector implementation of the ROL::Vector interface.
const Real getProbability(const int i) const
Provides the std::vector implementation of the ROL::Vector interface.
const ROL::Ptr< const ProbabilityVector< Real > > getProbabilityVector(void) const
const ROL::Ptr< const AtomVector< Real > > getAtomVector(void) const
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:84
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:167
constexpr auto dim