ROL
ROL_Triangle.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_TRIANGLE_HPP
45 #define ROL_TRIANGLE_HPP
46 
47 #include "ROL_Distribution.hpp"
48 #include "ROL_ParameterList.hpp"
49 
50 namespace ROL {
51 
52 template<class Real>
53 class Triangle : public Distribution<Real> {
54 private:
55  Real a_;
56  Real b_;
57  Real c_;
58 
59 public:
60  Triangle(const Real a = 0., const Real b = 0.5, const Real c = 1.)
61  : a_(std::min(a,std::min(b,c))),
62  b_(std::max(std::min(a,b),std::min(std::max(a,b),c))),
63  c_(std::max(a,std::max(b,c))) {}
64 
65  Triangle(ROL::ParameterList &parlist) {
66  Real a = parlist.sublist("SOL").sublist("Distribution").sublist("Triangle").get("Lower Bound",0.);
67  Real b = parlist.sublist("SOL").sublist("Distribution").sublist("Triangle").get("Peak Location",0.5);
68  Real c = parlist.sublist("SOL").sublist("Distribution").sublist("Triangle").get("Upper Bound",1.);
69  a_ = std::min(a,std::min(b,c));
70  b_ = std::max(std::min(a,b),std::min(std::max(a,b),c));
71  c_ = std::max(a,std::max(b,c));
72  }
73 
74  Real evaluatePDF(const Real input) const {
75  Real d1 = b_-a_, d2 = c_-b_, d = c_-a_;
76  return ((input >= a_ && input < b_) ? 2.0*(input-a_)/(d*d1) :
77  ((input >= b_ && input < c_) ? 2.0*(c_-input)/(d*d2) :
78  0.0));
79  }
80 
81  Real evaluateCDF(const Real input) const {
82  Real d1 = b_-a_, d2 = c_-b_, d = c_-a_;
83  return ((input < a_) ? 0.0 :
84  ((input >= a_ && input < b_) ?
85  std::pow(input-a_,2.0)/(d*d1) :
86  ((input >= b_ && input < c_) ?
87  1.0-std::pow(c_-input,2.0)/(d*d2) :
88  1.0)));
89  }
90 
91  Real integrateCDF(const Real input) const {
92  Real d1 = b_-a_, d2 = c_-b_, d = c_-a_;
93  return ((input < a_) ? 0.0 :
94  ((input >= a_ && input < b_) ?
95  std::pow(input-a_,3.0)/(3.0*d*d1) :
96  ((input >= b_ && input < c_) ?
97  d1*d1/(3.0*d)+(input-b_)+(std::pow(c_-input,3.0)-d2*d2*d2)/(3.0*d*d2) :
98  d1*d1/(3.0*d)+(input-b_)-d2*d2/(3.0*d))));
99  }
100 
101  Real invertCDF(const Real input) const {
102  Real d1 = b_-a_, d2 = c_-b_, d = c_-a_;
103  return ((input <= d1/d) ? a_ + std::sqrt(input*d1*d) :
104  c_ - std::sqrt((1.0-input)*d2*d));
105  }
106 
107  Real moment(const size_t m) const {
108  Real d1 = b_-a_, d2 = c_-b_, d = c_-a_;
109  Real am1 = std::pow(a_,m+1), am2 = a_*am1;
110  Real bm1 = std::pow(b_,m+1), bm2 = b_*bm1;
111  Real cm1 = std::pow(c_,m+1), cm2 = c_*cm1;
112  return (2./d)*(((bm2-am2)/((Real)m+2)-a_*(bm1-am1)/((Real)m+1))/d1
113  +(c_*(cm1-bm1)/((Real)m+1)-(cm2-bm2)/((Real)m+2))/d2);
114  }
115 
116  Real lowerBound(void) const {
117  return a_;
118  }
119 
120  Real upperBound(void) const {
121  return c_;
122  }
123 
124  void test(std::ostream &outStream = std::cout ) const {
125  size_t size = 7;
126  std::vector<Real> X(size,0.);
127  std::vector<int> T(size,0);
128  X[0] = a_-4.*(Real)rand()/(Real)RAND_MAX;
129  T[0] = 0;
130  X[1] = a_;
131  T[1] = 1;
132  X[2] = (b_-a_)*(Real)rand()/(Real)RAND_MAX + a_;
133  T[2] = 0;
134  X[3] = b_;
135  T[3] = 1;
136  X[4] = (c_-b_)*(Real)rand()/(Real)RAND_MAX + b_;
137  T[4] = 0;
138  X[5] = c_;
139  T[5] = 1;
140  X[6] = c_+4.*(Real)rand()/(Real)RAND_MAX;
141  T[6] = 0;
142  Distribution<Real>::test(X,T,outStream);
143  }
144 };
145 
146 }
147 
148 #endif
void test(std::ostream &outStream=std::cout) const
Triangle(const Real a=0., const Real b=0.5, const Real c=1.)
Real moment(const size_t m) const
Real upperBound(void) const
Real integrateCDF(const Real input) const
Real lowerBound(void) const
Real invertCDF(const Real input) const
virtual void test(std::ostream &outStream=std::cout) const
Real evaluatePDF(const Real input) const
Triangle(ROL::ParameterList &parlist)
Real evaluateCDF(const Real input) const