44 #ifndef ROL_TRUNCATEDGAUSSIAN_HPP 45 #define ROL_TRUNCATEDGAUSSIAN_HPP 49 #include "Teuchos_RCP.hpp" 50 #include "Teuchos_ParameterList.hpp" 62 Teuchos::RCP<Gaussian<Real> >
gauss_;
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) {
77 alpha_ = (a_-
mean_)/sdev_;
78 beta_ = (b_-
mean_)/sdev_;
79 phi_ = gauss_->evaluateCDF(alpha_);
80 Z_ = gauss_->evaluateCDF(beta_)-gauss_->evaluateCDF(alpha_);
84 const Real zero(0), one(1);
85 Teuchos::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);
91 b_ = std::max(b_,tmp);
93 mean_ = TGlist.get(
"Mean",zero);
94 sdev_ = TGlist.get(
"Standard Deviation",one);
95 sdev_ = (sdev_ > zero) ? sdev_ : one;
99 alpha_ = (a_-
mean_)/sdev_;
100 beta_ = (b_-
mean_)/sdev_;
101 phi_ = gauss_->evaluateCDF(alpha_);
102 Z_ = gauss_->evaluateCDF(beta_)-gauss_->evaluateCDF(alpha_);
106 const Real zero(0), xi = (input-
mean_)/sdev_;
107 return ((input <= a_) ? zero : ((input >=
b_) ? zero :
108 gauss_->evaluatePDF(xi)/(sdev_*
Z_)));
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_));
118 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
119 ">>> ERROR (ROL::TruncatedGaussian): Truncated Gaussian integrateCDF not implemented!");
124 const Real x = gauss_->invertCDF(Z_*input+phi_);
125 return sdev_*x +
mean_;
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_;
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;
139 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
140 ">>> ERROR (ROL::TruncatedGaussian): Truncated Gaussian moment not implemented for m > 2!");
153 void test(std::ostream &outStream = std::cout )
const {
155 std::vector<Real> X(size,0);
156 std::vector<int> T(size,0);
158 X[0] = a_-four*(Real)rand()/(Real)RAND_MAX;
162 X[2] = (b_-
a_)*(Real)rand()/(Real)RAND_MAX + a_;
166 X[4] = b_+four*(Real)rand()/(Real)RAND_MAX;
Real moment(const size_t m) const
Real invertCDF(const Real input) const
Teuchos::RCP< Gaussian< Real > > gauss_
void test(std::ostream &outStream=std::cout) const
virtual void test(std::ostream &outStream=std::cout) const
Real evaluateCDF(const Real input) const
Real integrateCDF(const Real input) const
Real lowerBound(void) const
TruncatedGaussian(Teuchos::ParameterList &parlist)
Real upperBound(void) const
TruncatedGaussian(const Real lo=-1., const Real up=1., const Real mean=0., const Real sdev=1.)
Real evaluatePDF(const Real input) const