ROL
ROL_Kumaraswamy.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_KUMARASWAMY_HPP
11 #define ROL_KUMARASWAMY_HPP
12 
13 #include "ROL_Distribution.hpp"
14 #include "ROL_ParameterList.hpp"
15 
16 #include <math.h>
17 
18 namespace ROL {
19 
20 template<class Real>
21 class Kumaraswamy : public Distribution<Real> {
22 private:
23  Real a_;
24  Real b_;
25  Real exp1_;
26  Real exp2_;
27  Real bGAMMAb_;
28 
29  size_t nchoosek(const size_t n, const size_t k) const {
30  return ((k==0) ? 1 : (n*nchoosek(n-1,k-1))/k);
31  }
32 
33 public:
34  Kumaraswamy(const Real a = 0., const Real b = 1.,
35  const Real exp1 = 0.5, const Real exp2 = 0.5)
36  : a_(std::min(a,b)), b_(std::max(a,b)),
37  exp1_((exp1>0.) ? exp1 : 0.5), exp2_((exp2>0.) ? exp2 : 0.5) {
38  bGAMMAb_ = exp2_*tgamma(exp2_);
39  }
40 
41  Kumaraswamy(ROL::ParameterList &parlist) {
42  a_ = parlist.sublist("SOL").sublist("Distribution").sublist("Kumaraswamy").get("Lower Bound",0.);
43  b_ = parlist.sublist("SOL").sublist("Distribution").sublist("Kumaraswamy").get("Upper Bound",1.);
44  Real tmp = a_;
45  a_ = std::min(a_,b_);
46  b_ = std::max(b_,tmp);
47 
48  exp1_ = parlist.sublist("SOL").sublist("Distribution").sublist("Kumaraswamy").get("Exponent 1",0.5);
49  exp2_ = parlist.sublist("SOL").sublist("Distribution").sublist("Kumaraswamy").get("Exponent 2",0.5);
50  exp1_ = (exp1_ > 0.) ? exp1_ : 0.5;
51  exp2_ = (exp2_ > 0.) ? exp2_ : 0.5;
52 
53  bGAMMAb_ = exp2_*tgamma(exp2_);
54  }
55 
56  Real evaluatePDF(const Real input) const {
57  Real x = (input - a_)/(b_-a_);
58  return ((x <= 0.) ? 0. : ((x >= 1.) ? 0. :
59  exp1_*exp2_*std::pow(x,exp1_-1)*std::pow(1.-std::pow(x,exp1_),exp2_-1)));
60  }
61 
62  Real evaluateCDF(const Real input) const {
63  Real x = (input - a_)/(b_-a_);
64  return ((x <= 0.) ? 0. : ((x >= 1.) ? 1. : 1.-std::pow(1.-std::pow(x,exp1_),exp2_)));
65  }
66 
67  Real integrateCDF(const Real input) const {
68  ROL_TEST_FOR_EXCEPTION( true, std::invalid_argument,
69  ">>> ERROR (ROL::Kumaraswamy): Kumaraswamy integrateCDF not implemented!");
70  }
71 
72  Real invertCDF(const Real input) const {
73  Real x = std::pow(1.-std::pow(1.-input,1./exp2_),1./exp1_);
74  return x*(b_-a_) + a_;
75  }
76 
77  Real moment(const size_t m) const {
78  Real val = 0., binom = 0., moment = 0.;
79  for (size_t i = 0; i < m+1; i++) {
80  moment = bGAMMAb_*tgamma(1.+(Real)i/exp1_)/tgamma(1.+exp2_+(Real)i/exp1_);
81  binom = (Real)nchoosek(m,i);
82  val += binom*std::pow(a_,m-i)*std::pow(b_-a_,i+1)*moment;
83  }
84  return val;
85  }
86 
87  Real lowerBound(void) const {
88  return a_;
89  }
90 
91  Real upperBound(void) const {
92  return b_;
93  }
94 
95  void test(std::ostream &outStream = std::cout ) const {
96  size_t size = 5;
97  std::vector<Real> X(size,0.);
98  std::vector<int> T(size,0);
99  X[0] = a_-4.*(Real)rand()/(Real)RAND_MAX;
100  T[0] = 0;
101  X[1] = a_;
102  T[1] = 1;
103  X[2] = (b_-a_)*(Real)rand()/(Real)RAND_MAX + a_;
104  T[2] = 0;
105  X[3] = b_;
106  T[3] = 1;
107  X[4] = b_+4.0*(Real)rand()/(Real)RAND_MAX;
108  T[4] = 0;
109  Distribution<Real>::test(X,T,outStream);
110  }
111 };
112 
113 }
114 
115 #endif
Kumaraswamy(ROL::ParameterList &parlist)
Kumaraswamy(const Real a=0., const Real b=1., const Real exp1=0.5, const Real exp2=0.5)
Real invertCDF(const Real input) const
Real integrateCDF(const Real input) const
Real moment(const size_t m) const
size_t nchoosek(const size_t n, const size_t k) const
Real lowerBound(void) const
virtual void test(std::ostream &outStream=std::cout) const
Real upperBound(void) const
Real evaluatePDF(const Real input) const
void test(std::ostream &outStream=std::cout) const
Real evaluateCDF(const Real input) const