44 #ifndef ROL_GAUSSIAN_HPP 45 #define ROL_GAUSSIAN_HPP 48 #include "Teuchos_ParameterList.hpp" 63 Real
erfi(
const Real p)
const {
64 const Real zero(0), half(0.5), one(1), two(2), pi(Teuchos::ScalarTraits<Real>::pi());
66 if ( std::abs(p) > static_cast<Real>(0.7) ) {
67 Real sgn = (p < zero) ? -one : one;
68 z = std::sqrt(-std::log((one-sgn*p)*half));
69 val = sgn*(((c_[3]*z+c_[2])*z+c_[1])*z + c_[0])/((d_[1]*z+d_[0])*z + one);
73 val = p*(((a_[3]*z+a_[2])*z+a_[1])*z + a_[0])/((((b_[3]*z+b_[2])*z+b_[1])*z+b_[0])*z+one);
75 val -= (erf(val)-p)/(two/std::sqrt(pi) * std::exp(-val*val));
76 val -= (erf(val)-p)/(two/std::sqrt(pi) * std::exp(-val*val));
82 Gaussian(
const Real mean = 0.,
const Real variance = 1.)
83 : mean_(mean), variance_((variance>0.) ? variance : 1.) {
84 a_.clear(); a_.resize(4,0.); b_.clear(); b_.resize(4,0.);
85 c_.clear(); c_.resize(4,0.); d_.clear(); d_.resize(2,0.);
86 a_[0] = 0.886226899; a_[1] = -1.645349621; a_[2] = 0.914624893; a_[3] = -0.140543331;
87 b_[0] = -2.118377725; b_[1] = 1.442710462; b_[2] = -0.329097515; b_[3] = 0.012229801;
88 c_[0] = -1.970840454; c_[1] = -1.624906493; c_[2] = 3.429567803; c_[3] = 1.641345311;
89 d_[0] = 3.543889200; d_[1] = 1.637067800;
93 mean_ = parlist.sublist(
"SOL").sublist(
"Distribution").sublist(
"Gaussian").get(
"Mean",0.);
94 variance_ = parlist.sublist(
"SOL").sublist(
"Distribution").sublist(
"Gaussian").get(
"Variance",1.);
95 variance_ = (variance_ > 0.) ?
variance_ : 1.;
96 a_.clear(); a_.resize(4,0.); b_.clear(); b_.resize(4,0.);
97 c_.clear(); c_.resize(4,0.); d_.clear(); d_.resize(2,0.);
98 a_[0] = 0.886226899; a_[1] = -1.645349621; a_[2] = 0.914624893; a_[3] = -0.140543331;
99 b_[0] = -2.118377725; b_[1] = 1.442710462; b_[2] = -0.329097515; b_[3] = 0.012229801;
100 c_[0] = -1.970840454; c_[1] = -1.624906493; c_[2] = 3.429567803; c_[3] = 1.641345311;
101 d_[0] = 3.543889200; d_[1] = 1.637067800;
105 return std::exp(-std::pow(input-mean_,2)/(2.*variance_))/(std::sqrt(2.*Teuchos::ScalarTraits<Real>::pi()*variance_));
109 const Real half(0.5), one(1), two(2);
110 return half*(one+erf((input-mean_)/std::sqrt(two*variance_)));
114 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
115 ">>> ERROR (ROL::Gaussian): Gaussian integrateCDF not implemented!");
116 return ((input < mean_) ? 0.0 : input);
121 const Real zero(0), half(0.5), one(1), eps(ROL_EPSILON<Real>());
122 if ( input <= eps ) {
125 if ( input >= one-eps ) {
128 Real a = eps, b = one-eps, c = zero;
131 Real sa = ((fa < zero) ? -one : ((fa > zero) ? one : zero));
133 for (
size_t i = 0; i < 100; i++) {
136 sc = ((fc < zero) ? -one : ((fc > zero) ? one : zero));
137 if ( fc == zero || (b-a)*half < eps ) {
140 if ( sc == sa ) { a = c; fa = fc; sa = sc; }
149 case 1: val =
mean_;
break;
150 case 2: val = std::pow(mean_,2) +
variance_;
break;
151 case 3: val = std::pow(mean_,3)
153 case 4: val = std::pow(mean_,4)
154 + 6.*std::pow(mean_,2)*variance_
155 + 3.*std::pow(variance_,2);
break;
156 case 5: val = std::pow(mean_,5)
157 + 10.*std::pow(mean_,3)*variance_
158 + 15.*mean_*std::pow(variance_,2);
break;
159 case 6: val = std::pow(mean_,6)
160 + 15.*std::pow(mean_,4)*variance_
161 + 45.*std::pow(mean_*variance_,2)
162 + 15.*std::pow(variance_,3);
break;
163 case 7: val = std::pow(mean_,7)
164 + 21.*std::pow(mean_,5)*variance_
165 + 105.*std::pow(mean_,3)*std::pow(variance_,2)
166 + 105.*mean_*std::pow(variance_,3);
break;
167 case 8: val = std::pow(mean_,8)
168 + 28.*std::pow(mean_,6)*variance_
169 + 210.*std::pow(mean_,4)*std::pow(variance_,2)
170 + 420.*std::pow(mean_,2)*std::pow(variance_,3)
171 + 105.*std::pow(variance_,4);
break;
173 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
174 ">>> ERROR (ROL::Distribution): Gaussian moment not implemented for m > 8!");
180 return ROL_NINF<Real>();
184 return ROL_INF<Real>();
187 void test(std::ostream &outStream = std::cout )
const {
189 std::vector<Real> X(size,4.*(Real)rand()/(Real)RAND_MAX - 2.);
190 std::vector<int> T(size,0);
Real evaluateCDF(const Real input) const
Real evaluatePDF(const Real input) const
virtual void test(std::ostream &outStream=std::cout) const
Real upperBound(void) const
void test(std::ostream &outStream=std::cout) const
Gaussian(Teuchos::ParameterList &parlist)
Real moment(const size_t m) const
Gaussian(const Real mean=0., const Real variance=1.)
Real lowerBound(void) const
Real invertCDF(const Real input) const
Real integrateCDF(const Real input) const
Real erfi(const Real p) const