ROL
ROL_PointwiseCDFObjective.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_POINTWISECDFOBJECTIVE_H
45 #define ROL_POINTWISECDFOBJECTIVE_H
46 
47 #include "ROL_Objective.hpp"
48 #include "ROL_BatchManager.hpp"
49 #include "ROL_Vector.hpp"
50 #include "ROL_Distribution.hpp"
51 #include "ROL_Ptr.hpp"
52 #include <math.h>
53 
54 namespace ROL {
55 
56 template <class Real>
57 class PointwiseCDFObjective : public Objective<Real> {
58 private:
59  std::vector<ROL::Ptr<Distribution<Real> > > dist_;
60  ROL::Ptr<BatchManager<Real> > bman_;
61  const Real scale_;
62  const Real sqrt2_;
63  const Real sqrtpi_;
64 
65  Real valueCDF(const int dim, const Real loc, const SROMVector<Real> &x) const {
66  const int numSamples = x.getNumSamples();
67  Real val = 0., hs = 0., xpt = 0., xwt = 0.;
68  for (int k = 0; k < numSamples; k++) {
69  xpt = (*x.getPoint(k))[dim]; xwt = x.getWeight(k);
70  hs = 0.5 * (1. + erf((loc-xpt)/(sqrt2_*scale_)));
71  val += xwt * hs;
72  }
73  return val;
74  }
75 
76  Real gradientCDF(std::vector<Real> &gradx, std::vector<Real> &gradp,
77  const int dim, const Real loc, const SROMVector<Real> &x) const {
78  const int numSamples = x.getNumSamples();
79  gradx.resize(numSamples,0.); gradp.resize(numSamples,0.);
80  Real val = 0., hs = 0., xpt = 0., xwt = 0.;
81  for (int k = 0; k < numSamples; k++) {
82  xpt = (*x.getPoint(k))[dim]; xwt = x.getWeight(k);
83  hs = 0.5 * (1. + erf((loc-xpt)/(sqrt2_*scale_)));
84  val += xwt * hs;
85  gradx[k] = -(xwt/(sqrt2_*sqrtpi_*scale_))
86  * std::exp(-std::pow((loc-xpt)/(sqrt2_*scale_),2));
87  gradp[k] = hs;
88  }
89  return val;
90  }
91 
92  Real hessVecCDF(std::vector<Real> &hvx,
93  const int dim, const Real loc, const SROMVector<Real> &x, const SROMVector<Real> &v) const {
94  const int numSamples = x.getNumSamples();
95  hvx.resize(numSamples,0.);
96  Real val = 0., hs = 0., xpt = 0., xwt = 0., scale3 = std::pow(scale_,3);
97  for (int k = 0; k < numSamples; k++) {
98  xpt = (*x.getPoint(k))[dim]; xwt = x.getWeight(k);
99  hs = 0.5 * (1. + erf((loc-xpt)/(sqrt2_*scale_)));
100  val += xwt * hs;
101  hvx[k] = -(xwt/(sqrt2_*sqrtpi_*scale3))
102  * std::exp(-std::pow((loc-xpt)/(sqrt2_*scale_),2)) * (loc-xpt);
103  }
104  return val;
105  }
106 
107 public:
108  PointwiseCDFObjective(const std::vector<ROL::Ptr<Distribution<Real> > > &dist,
109  ROL::Ptr<BatchManager<Real> > &bman,
110  const Real scale = 1.e-2)
111  : Objective<Real>(), dist_(dist), bman_(bman), scale_(scale),
112  sqrt2_(std::sqrt(2.)), sqrtpi_(std::sqrt(ROL::ScalarTraits<Real>::pi())) {}
113 
114  Real value( const Vector<Real> &x, Real &tol ) {
115  const SROMVector<Real> &ex = dynamic_cast<const SROMVector<Real>&>(x);
116  const int dimension = ex.getDimension();
117  const int numSamples = ex.getNumSamples();
118  Real val = 0., diff = 0., xpt = 0., sum = 0.;
119  for (int d = 0; d < dimension; d++) {
120  for (int k = 0; k < numSamples; k++) {
121  xpt = (*ex.getPoint(k))[d];
122  diff = (valueCDF(d,xpt,ex)-dist_[d]->evaluateCDF(xpt));
123  val += std::pow(diff,2);
124  }
125  }
126  bman_->sumAll(&val,&sum,1);
127  return 0.5*sum;
128  }
129 
130  void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
131  SROMVector<Real> &eg = dynamic_cast<SROMVector<Real>&>(g);
132  const SROMVector<Real> &ex = dynamic_cast<const SROMVector<Real>&>(x);
133  const int dimension = ex.getDimension();
134  const int numSamples = ex.getNumSamples();
135  std::vector<Real> gradx(numSamples,0.), gradp(numSamples,0.);
136  Real diff = 0., xpt = 0., val = 0., sum = 0.;
137  std::vector<Real> val_wt(numSamples,0.), tmp(dimension,0.);
138  std::vector<std::vector<Real> > val_pt(numSamples,tmp);
139  for (int d = 0; d < dimension; d++) {
140  for (int k = 0; k < numSamples; k++) {
141  xpt = (*ex.getPoint(k))[d];
142  val = gradientCDF(gradx,gradp,d,xpt,ex);
143  diff = (val-dist_[d]->evaluateCDF(xpt));
144  sum = 0.;
145  for (int j = 0; j < numSamples; j++) {
146  (val_pt[j])[d] += diff * gradx[j];
147  val_wt[j] += diff * gradp[j];
148  sum -= gradx[j];
149  }
150  (val_pt[k])[d] += diff * (sum - dist_[d]->evaluatePDF(xpt));
151  }
152  }
153  for (int k = 0; k < numSamples; k++) {
154  eg.setPoint(k,val_pt[k]);
155  eg.setWeight(k,val_wt[k]);
156  }
157  }
158 
159 // void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
160 // }
161 }; // class LinearCombinationObjective
162 
163 } // namespace ROL
164 
165 #endif
Provides the interface to evaluate objective functions.
Real valueCDF(const int dim, const Real loc, const SROMVector< Real > &x) const
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
PointwiseCDFObjective(const std::vector< ROL::Ptr< Distribution< Real > > > &dist, ROL::Ptr< BatchManager< Real > > &bman, const Real scale=1.e-2)
Real value(const Vector< Real > &x, Real &tol)
Compute value.
Real hessVecCDF(std::vector< Real > &hvx, const int dim, const Real loc, const SROMVector< Real > &x, const SROMVector< Real > &v) const
Provides the std::vector implementation of the ROL::Vector interface.
ROL::Ptr< BatchManager< Real > > bman_
Real gradientCDF(std::vector< Real > &gradx, std::vector< Real > &gradp, const int dim, const Real loc, const SROMVector< Real > &x) const
constexpr auto dim
std::vector< ROL::Ptr< Distribution< Real > > > dist_