ROL
ROL_TruncatedGaussian.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_TRUNCATEDGAUSSIAN_HPP
11 #define ROL_TRUNCATEDGAUSSIAN_HPP
12 
13 #include "ROL_Distribution.hpp"
14 #include "ROL_Gaussian.hpp"
15 #include "ROL_Ptr.hpp"
16 #include "ROL_ParameterList.hpp"
17 
18 namespace ROL {
19 
20 template<class Real>
21 class TruncatedGaussian : public Distribution<Real> {
22 private:
23  Real a_;
24  Real b_;
25  Real mean_;
26  Real sdev_;
27 
28  ROL::Ptr<Gaussian<Real> > gauss_;
29 
30  Real alpha_;
31  Real beta_;
32  Real phi_;
33  Real Z_;
34 
35 public:
36 
37  TruncatedGaussian(const Real lo = -1., const Real up = 1.,
38  const Real mean = 0., const Real sdev = 1.)
39  : a_((lo < up) ? lo : up), b_((up > lo) ? up : lo),
40  mean_(mean), sdev_((sdev>0) ? sdev : 1) {
41  //Real var = sdev_*sdev_;
42  gauss_ = ROL::makePtr<Gaussian<Real>>();
43  alpha_ = (a_-mean_)/sdev_;
44  beta_ = (b_-mean_)/sdev_;
45  phi_ = gauss_->evaluateCDF(alpha_);
46  Z_ = gauss_->evaluateCDF(beta_)-gauss_->evaluateCDF(alpha_);
47  }
48 
49  TruncatedGaussian(ROL::ParameterList &parlist) {
50  const Real zero(0), one(1);
51  ROL::ParameterList TGlist
52  = parlist.sublist("SOL").sublist("Distribution").sublist("Truncated Gaussian");
53  a_ = TGlist.get("Lower Bound",-one);
54  b_ = TGlist.get("Upper Bound", one);
55  Real tmp = a_;
56  a_ = std::min(a_,b_);
57  b_ = std::max(b_,tmp);
58 
59  mean_ = TGlist.get("Mean",zero);
60  sdev_ = TGlist.get("Standard Deviation",one);
61  sdev_ = (sdev_ > zero) ? sdev_ : one;
62 
63  //Real var = sdev_*sdev_;
64  gauss_ = ROL::makePtr<Gaussian<Real>>();
65  alpha_ = (a_-mean_)/sdev_;
66  beta_ = (b_-mean_)/sdev_;
67  phi_ = gauss_->evaluateCDF(alpha_);
68  Z_ = gauss_->evaluateCDF(beta_)-gauss_->evaluateCDF(alpha_);
69  }
70 
71  Real evaluatePDF(const Real input) const {
72  const Real zero(0), xi = (input-mean_)/sdev_;
73  return ((input <= a_) ? zero : ((input >= b_) ? zero :
74  gauss_->evaluatePDF(xi)/(sdev_*Z_)));
75  }
76 
77  Real evaluateCDF(const Real input) const {
78  const Real zero(0), one(1), xi = (input-mean_)/sdev_;
79  return ((input <= a_) ? zero : ((input >= b_) ? one :
80  (gauss_->evaluateCDF(xi)-phi_)/Z_));
81  }
82 
83  Real integrateCDF(const Real input) const {
84  ROL_TEST_FOR_EXCEPTION( true, std::invalid_argument,
85  ">>> ERROR (ROL::TruncatedGaussian): Truncated Gaussian integrateCDF not implemented!");
86  //return ((input < 0.5*(a_+b_)) ? 0.0 : input - 0.5*(a_+b_));
87  }
88 
89  Real invertCDF(const Real input) const {
90  const Real x = gauss_->invertCDF(Z_*input+phi_);
91  return sdev_*x + mean_;
92  }
93 
94  Real moment(const size_t m) const {
95  const Real phiA = gauss_->evaluatePDF(alpha_);
96  const Real phiB = gauss_->evaluatePDF(beta_);
97  const Real mean = mean_ + sdev_*(phiA-phiB)/Z_;
98  const Real var = sdev_*sdev_;
99  const Real one(1);
100  Real val(0);
101  switch(m) {
102  case 1: val = mean; break;
103  case 2: val = var*(one+(alpha_*phiA-beta_*phiB)/Z_-std::pow((phiA-phiB)/Z_,2))+mean*mean; break;
104  default:
105  ROL_TEST_FOR_EXCEPTION( true, std::invalid_argument,
106  ">>> ERROR (ROL::TruncatedGaussian): Truncated Gaussian moment not implemented for m > 2!");
107  }
108  return val;
109  }
110 
111  Real lowerBound(void) const {
112  return a_;
113  }
114 
115  Real upperBound(void) const {
116  return b_;
117  }
118 
119  void test(std::ostream &outStream = std::cout ) const {
120  size_t size = 5;
121  std::vector<Real> X(size,0);
122  std::vector<int> T(size,0);
123  const Real four(4);
124  X[0] = a_-four*(Real)rand()/(Real)RAND_MAX;
125  T[0] = 0;
126  X[1] = a_;
127  T[1] = 1;
128  X[2] = (b_-a_)*(Real)rand()/(Real)RAND_MAX + a_;
129  T[2] = 0;
130  X[3] = b_;
131  T[3] = 1;
132  X[4] = b_+four*(Real)rand()/(Real)RAND_MAX;
133  T[4] = 0;
134  Distribution<Real>::test(X,T,outStream);
135  }
136 };
137 
138 }
139 
140 #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.)