ROL
ROL_TruncatedGaussian.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_TRUNCATEDGAUSSIAN_HPP
45 #define ROL_TRUNCATEDGAUSSIAN_HPP
46 
47 #include "ROL_Distribution.hpp"
48 #include "ROL_Gaussian.hpp"
49 #include "ROL_Ptr.hpp"
50 #include "ROL_ParameterList.hpp"
51 
52 namespace ROL {
53 
54 template<class Real>
55 class TruncatedGaussian : public Distribution<Real> {
56 private:
57  Real a_;
58  Real b_;
59  Real mean_;
60  Real sdev_;
61 
62  ROL::Ptr<Gaussian<Real> > gauss_;
63 
64  Real alpha_;
65  Real beta_;
66  Real phi_;
67  Real Z_;
68 
69 public:
70 
71  TruncatedGaussian(const Real lo = -1., const Real up = 1.,
72  const Real mean = 0., const Real sdev = 1.)
73  : a_((lo < up) ? lo : up), b_((up > lo) ? up : lo),
74  mean_(mean), sdev_((sdev>0) ? sdev : 1) {
75  //Real var = sdev_*sdev_;
76  gauss_ = ROL::makePtr<Gaussian<Real>>();
77  alpha_ = (a_-mean_)/sdev_;
78  beta_ = (b_-mean_)/sdev_;
79  phi_ = gauss_->evaluateCDF(alpha_);
80  Z_ = gauss_->evaluateCDF(beta_)-gauss_->evaluateCDF(alpha_);
81  }
82 
83  TruncatedGaussian(ROL::ParameterList &parlist) {
84  const Real zero(0), one(1);
85  ROL::ParameterList TGlist
86  = parlist.sublist("SOL").sublist("Distribution").sublist("Truncated Gaussian");
87  a_ = TGlist.get("Lower Bound",-one);
88  b_ = TGlist.get("Upper Bound", one);
89  Real tmp = a_;
90  a_ = std::min(a_,b_);
91  b_ = std::max(b_,tmp);
92 
93  mean_ = TGlist.get("Mean",zero);
94  sdev_ = TGlist.get("Standard Deviation",one);
95  sdev_ = (sdev_ > zero) ? sdev_ : one;
96 
97  //Real var = sdev_*sdev_;
98  gauss_ = ROL::makePtr<Gaussian<Real>>();
99  alpha_ = (a_-mean_)/sdev_;
100  beta_ = (b_-mean_)/sdev_;
101  phi_ = gauss_->evaluateCDF(alpha_);
102  Z_ = gauss_->evaluateCDF(beta_)-gauss_->evaluateCDF(alpha_);
103  }
104 
105  Real evaluatePDF(const Real input) const {
106  const Real zero(0), xi = (input-mean_)/sdev_;
107  return ((input <= a_) ? zero : ((input >= b_) ? zero :
108  gauss_->evaluatePDF(xi)/(sdev_*Z_)));
109  }
110 
111  Real evaluateCDF(const Real input) const {
112  const Real zero(0), one(1), xi = (input-mean_)/sdev_;
113  return ((input <= a_) ? zero : ((input >= b_) ? one :
114  (gauss_->evaluateCDF(xi)-phi_)/Z_));
115  }
116 
117  Real integrateCDF(const Real input) const {
118  ROL_TEST_FOR_EXCEPTION( true, std::invalid_argument,
119  ">>> ERROR (ROL::TruncatedGaussian): Truncated Gaussian integrateCDF not implemented!");
120  //return ((input < 0.5*(a_+b_)) ? 0.0 : input - 0.5*(a_+b_));
121  }
122 
123  Real invertCDF(const Real input) const {
124  const Real x = gauss_->invertCDF(Z_*input+phi_);
125  return sdev_*x + mean_;
126  }
127 
128  Real moment(const size_t m) const {
129  const Real phiA = gauss_->evaluatePDF(alpha_);
130  const Real phiB = gauss_->evaluatePDF(beta_);
131  const Real mean = mean_ + sdev_*(phiA-phiB)/Z_;
132  const Real var = sdev_*sdev_;
133  const Real one(1);
134  Real val(0);
135  switch(m) {
136  case 1: val = mean; break;
137  case 2: val = var*(one+(alpha_*phiA-beta_*phiB)/Z_-std::pow((phiA-phiB)/Z_,2))+mean*mean; break;
138  default:
139  ROL_TEST_FOR_EXCEPTION( true, std::invalid_argument,
140  ">>> ERROR (ROL::TruncatedGaussian): Truncated Gaussian moment not implemented for m > 2!");
141  }
142  return val;
143  }
144 
145  Real lowerBound(void) const {
146  return a_;
147  }
148 
149  Real upperBound(void) const {
150  return b_;
151  }
152 
153  void test(std::ostream &outStream = std::cout ) const {
154  size_t size = 5;
155  std::vector<Real> X(size,0);
156  std::vector<int> T(size,0);
157  const Real four(4);
158  X[0] = a_-four*(Real)rand()/(Real)RAND_MAX;
159  T[0] = 0;
160  X[1] = a_;
161  T[1] = 1;
162  X[2] = (b_-a_)*(Real)rand()/(Real)RAND_MAX + a_;
163  T[2] = 0;
164  X[3] = b_;
165  T[3] = 1;
166  X[4] = b_+four*(Real)rand()/(Real)RAND_MAX;
167  T[4] = 0;
168  Distribution<Real>::test(X,T,outStream);
169  }
170 };
171 
172 }
173 
174 #endif
void test(std::ostream &outStream=std::cout) const
Real integrateCDF(const Real input) const
Real invertCDF(const Real input) const
Real evaluateCDF(const Real input) const
TruncatedGaussian(ROL::ParameterList &parlist)
Real moment(const size_t m) const
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Real evaluatePDF(const Real input) const
ROL::Ptr< Gaussian< Real > > gauss_
virtual void test(std::ostream &outStream=std::cout) const
TruncatedGaussian(const Real lo=-1., const Real up=1., const Real mean=0., const Real sdev=1.)