44 #ifndef ROL_SMOOTHEDPOE_HPP
45 #define ROL_SMOOTHEDPOE_HPP
86 const Real one(1), two(2);
89 Real ex = std::exp(-two*x/
eps_);
92 else if (deriv == 1) {
93 Real ex = std::exp(-two*x/
eps_);
94 val = (two/
eps_)*ex/std::pow(one+ex,2);
96 else if (deriv == 2) {
97 Real ex = std::exp(two*x/
eps_);
98 val = std::pow(two/
eps_,2)*ex*(one-ex)/std::pow(one+ex,3);
110 ROL::ParameterList &list = parlist.sublist(
"SOL").sublist(
"Probability").sublist(
"Smoothed POE");
112 eps_ = list.get<Real>(
"Smoothing Parameter");
117 const std::vector<Real> &xstat,
121 if ( std::abs(sp) > ROL_EPSILON<Real>() ) {
127 const std::vector<Real> &xstat,
136 const std::vector<Real> &xstat,
140 if ( std::abs(sp) > ROL_EPSILON<Real>() ) {
147 std::vector<Real> &gstat,
149 const std::vector<Real> &xstat,
156 const std::vector<Real> &vstat,
158 const std::vector<Real> &xstat,
163 if ( std::abs(sp1) > ROL_EPSILON<Real>() ) {
168 if ( std::abs(sp2) > ROL_EPSILON<Real>() ) {
176 std::vector<Real> &hvstat,
178 const std::vector<Real> &vstat,
180 const std::vector<Real> &xstat,
Provides the interface to evaluate objective functions.
void computeHessVec(Vector< Real > &hv, Objective< Real > &obj, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Real computeValue(Objective< Real > &obj, const Vector< Real > &x, Real &tol)
void getGradient(Vector< Real > &g, std::vector< Real > &gstat, const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure (sub)gradient.
Ptr< Vector< Real > > hv_
void updateValue(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal storage for value computation.
void updateHessVec(Objective< Real > &obj, const Vector< Real > &v, const std::vector< Real > &vstat, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal risk measure storage for Hessian-time-a-vector computation.
Ptr< Vector< Real > > dualVector_
Defines the linear algebra or vector space interface.
void sumAll(Real *input, Real *output, int dim) const
SmoothedPOE(ROL::ParameterList &parlist)
Real getValue(const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure value.
void getHessVec(Vector< Real > &hv, std::vector< Real > &hvstat, const Vector< Real > &v, const std::vector< Real > &vstat, const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure Hessian-times-a-vector.
void computeGradient(Vector< Real > &g, Objective< Real > &obj, const Vector< Real > &x, Real &tol)
Real smoothHeaviside(const Real x, const int deriv=0) const
SmoothedPOE(const Real threshold, const Real eps)
Real computeGradVec(Vector< Real > &g, Objective< Real > &obj, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Provides the interface to implement any functional that maps a random variable to a (extended) real n...
void updateGradient(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal risk measure storage for gradient computation.
Provides the implementation of the smoothed probability of exceedance.