ROL
ROL_Exponential.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_EXPONENTIAL_HPP
45 #define ROL_EXPONENTIAL_HPP
46 
47 #include "ROL_Distribution.hpp"
48 #include "ROL_ParameterList.hpp"
49 
50 namespace ROL {
51 
52 template<class Real>
53 class Exponential : public Distribution<Real> {
54 private:
55  Real loc_;
56  Real scale_;
57 
58  size_t compute_coeff(const size_t m, const size_t k) const {
59  if ( k == 0 || m == 0 || m == 1 ) {
60  return 1;
61  }
62  size_t val = 1;
63  for (size_t i = m-k; i < m; i++) {
64  val *= (i+1);
65  }
66  return val;
67  }
68 
69 public:
70  Exponential(const Real loc = 0., const Real scale = 1.)
71  : loc_(loc), scale_((scale>0.) ? scale : 1.) {}
72 
73  Exponential(ROL::ParameterList &parlist) {
74  loc_ = parlist.sublist("SOL").sublist("Distribution").sublist("Exponential").get("Location",0.);
75  scale_ = parlist.sublist("SOL").sublist("Distribution").sublist("Exponential").get("Scale",1.);
76  scale_ = (scale_ > 0.) ? scale_ : 1.;
77  }
78 
79  Real evaluatePDF(const Real input) const {
80  return ((input >= loc_) ? scale_*std::exp(-scale_*(input-loc_)) : 0.);
81  }
82 
83  Real evaluateCDF(const Real input) const {
84  return ((input >= loc_) ? 1.-std::exp(-scale_*(input-loc_)) : 0.);
85  }
86 
87  Real integrateCDF(const Real input) const {
88  return ((input >= loc_) ?
89  (input-loc_) - (1.-std::exp(-scale_*(input-loc_)))/scale_ : 0.);
90  }
91 
92  Real invertCDF(const Real input) const {
93  return -std::log(1.-input)/scale_;
94  }
95 
96  Real moment(const size_t m) const {
97  Real val = 0., coeff = 0.;
98  for (size_t i = 0; i < m+1; i++) {
99  coeff = compute_coeff(m,i);
100  val += coeff*std::pow(loc_,(Real)(m-i))/std::pow(scale_,(Real)i);
101  }
102  return val;
103  }
104 
105  Real lowerBound(void) const {
106  return 0.;
107  }
108 
109  Real upperBound(void) const {
110  return ROL_INF<Real>();
111  }
112 
113  void test(std::ostream &outStream = std::cout ) const {
114  size_t size = 3;
115  std::vector<Real> X(size,0.);
116  std::vector<int> T(size,0);
117  X[0] = loc_-4.0*(Real)rand()/(Real)RAND_MAX;
118  T[0] = 0;
119  X[1] = loc_;
120  T[1] = 1;
121  X[2] = loc_+4.0*(Real)rand()/(Real)RAND_MAX;
122  T[2] = 0;
123  Distribution<Real>::test(X,T,outStream);
124  }
125 };
126 
127 }
128 
129 #endif
Real integrateCDF(const Real input) const
Exponential(ROL::ParameterList &parlist)
Real invertCDF(const Real input) const
Real upperBound(void) const
Real evaluatePDF(const Real input) const
Real evaluateCDF(const Real input) const
size_t compute_coeff(const size_t m, const size_t k) const
Real moment(const size_t m) const
virtual void test(std::ostream &outStream=std::cout) const
Exponential(const Real loc=0., const Real scale=1.)
Real lowerBound(void) const
void test(std::ostream &outStream=std::cout) const