ROL
ROL_Smale.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_SMALE_HPP
45#define ROL_SMALE_HPP
46
47#include "ROL_Distribution.hpp"
48#include "ROL_ParameterList.hpp"
49
50namespace ROL {
51
52template<class Real>
53class Smale : public Distribution<Real> {
54private:
55 Real a_;
56 Real b_;
57
58public:
59 Smale(const Real a = 0., const Real b = 1.)
60 : a_(std::min(a,b)), b_(std::max(a,b)) {}
61
62 Smale(ROL::ParameterList &parlist) {
63 a_ = parlist.sublist("SOL").sublist("Distribution").sublist("Smale").get("Lower Bound",0.);
64 b_ = parlist.sublist("SOL").sublist("Distribution").sublist("Smale").get("Upper Bound",1.);
65 Real tmp = a_;
66 a_ = std::min(a_,b_);
67 b_ = std::max(b_,tmp);
68 }
69
70 Real evaluatePDF(const Real input) const {
71 Real val = std::pow(input-a_,2)+4.*b_*b_;
72 Real root = std::sqrt(val);
73 return 2.0*b_*b_/(val*root);
74 }
75
76 Real evaluateCDF(const Real input) const {
77 Real val = std::pow(input-a_,2)+4.*b_*b_;
78 Real root = std::sqrt(val);
79 return 0.5*(1.0+input/root);
80 }
81
82 Real integrateCDF(const Real input) const {
83 Real val = std::pow(input-a_,2)+4.*b_*b_;
84 Real root = std::sqrt(val);
85 return 0.5*(input+root);
86 }
87
88 Real invertCDF(const Real input) const {
89 Real x = a_;
90 Real fx = evaluateCDF(x)-input;
91 Real s = 0.0;
92 Real xs = 0.0;
93 Real a = 1.0;
94 Real tmp = 0.0;
95 for (int i = 0; i < 100; i++) {
96 if ( std::abs(fx) < ROL_EPSILON<Real>() ) { break; }
97 s = -fx/evaluatePDF(x);
98 a = 1.0;
99 xs = x + a*s;
100 tmp = fx;
101 fx = evaluateCDF(xs)-input;
102 while ( std::abs(fx) > (1.0 - 1.e-4*a)*std::abs(tmp) ) {
103 a *= 0.5;
104 xs = x + a*s;
105 fx = evaluateCDF(xs)-input;
106 }
107 x = xs;
108 }
109 return x;
110 }
111
112 Real moment(const size_t m) const {
113 ROL_TEST_FOR_EXCEPTION( true, std::invalid_argument,
114 ">>> ERROR (ROL::Smale): Smale moment is not implemented!");
115 }
116
117 Real lowerBound(void) const {
118 return ROL_NINF<Real>();
119 }
120
121 Real upperBound(void) const {
122 return ROL_INF<Real>();
123 }
124
125 void test(std::ostream &outStream = std::cout ) const {
126 size_t size = 1;
127 std::vector<Real> X(size,4.*(Real)rand()/(Real)RAND_MAX - 2.);
128 std::vector<int> T(size,0);
129 Distribution<Real>::test(X,T,outStream);
130 }
131};
132
133}
134
135#endif
virtual void test(std::ostream &outStream=std::cout) const
Real upperBound(void) const
Real lowerBound(void) const
Smale(ROL::ParameterList &parlist)
Definition ROL_Smale.hpp:62
Smale(const Real a=0., const Real b=1.)
Definition ROL_Smale.hpp:59
void test(std::ostream &outStream=std::cout) const
Real integrateCDF(const Real input) const
Definition ROL_Smale.hpp:82
Real evaluateCDF(const Real input) const
Definition ROL_Smale.hpp:76
Real evaluatePDF(const Real input) const
Definition ROL_Smale.hpp:70
Real moment(const size_t m) const
Real invertCDF(const Real input) const
Definition ROL_Smale.hpp:88