10 #ifndef ROL_MEANDEVIATIONFROMTARGET_HPP
11 #define ROL_MEANDEVIATIONFROMTARGET_HPP
14 #include "ROL_ParameterList.hpp"
54 std::vector<Ptr<Vector<Real> > >
pg0_;
55 std::vector<Ptr<Vector<Real> > >
pg_;
56 std::vector<Ptr<Vector<Real> > >
phv_;
86 ROL_TEST_FOR_EXCEPTION((oSize!=cSize),std::invalid_argument,
87 ">>> ERROR (ROL::MeanDeviationFromTarget): Order and coefficient arrays have different sizes!");
88 ROL_TEST_FOR_EXCEPTION((oSize!=tSize),std::invalid_argument,
89 ">>> ERROR (ROL::MeanDeviationFromTarget): Order and target arrays have different sizes!");
91 for (
int i = 0; i < oSize; i++) {
92 ROL_TEST_FOR_EXCEPTION((
order_[i] < two), std::invalid_argument,
93 ">>> ERROR (ROL::MeanDeviationFromTarget): Element of order array out of range!");
94 ROL_TEST_FOR_EXCEPTION((
coeff_[i] <
zero), std::invalid_argument,
95 ">>> ERROR (ROL::MeanDeviationFromTarget): Element of coefficient array out of range!");
98 ">>> ERROR (ROL::MeanDeviationFromTarget): PositiveFunction pointer is null!");
134 const std::vector<Real> &order,
135 const std::vector<Real> &coeff,
139 for (
uint i = 0; i < target.size(); i++ ) {
142 for (
uint i = 0; i < order.size(); i++ ) {
143 order_.push_back(order[i]);
145 for (
uint i = 0; i < coeff.size(); i++ ) {
146 coeff_.push_back(coeff[i]);
166 ROL::ParameterList &list
167 = parlist.sublist(
"SOL").sublist(
"Risk Measure").sublist(
"Mean Plus Deviation From Target");
170 target_ = ROL::getArrayFromStringParameter<double>(list,
"Targets");
172 order_ = ROL::getArrayFromStringParameter<double>(list,
"Orders");
174 coeff_ = ROL::getArrayFromStringParameter<double>(list,
"Coefficients");
177 std::string type = list.get<std::string>(
"Deviation Type");
178 if ( type ==
"Upper" ) {
181 else if ( type ==
"Absolute" ) {
185 ROL_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
186 ">>> (ROL::MeanDeviation): Deviation type is not recoginized!");
212 const std::vector<Real> &xstat,
214 Real diff(0), pf0(0);
225 const std::vector<Real> &xstat,
240 const std::vector<Real> &xstat,
242 Real diff(0), pf0(0), pf1(0), c(0), one(1);
249 c = std::pow(pf0,
order_[p]-one) * pf1;
257 std::vector<Real> &gstat,
259 const std::vector<Real> &xstat,
261 const Real
zero(0), one(1);
266 if ( pval_sum[p] >
zero ) {
276 const std::vector<Real> &vstat,
278 const std::vector<Real> &xstat,
280 const Real one(1), two(2);
281 Real diff(0), pf0(0), pf1(0), pf2(0), p0(0), p1(0), p2(0), c(0);
290 p0 = std::pow(pf0,
order_[p]);
291 p1 = std::pow(pf0,
order_[p]-one);
292 p2 = std::pow(pf0,
order_[p]-two);
293 c = -(
order_[p]-one)*p1*pf1;
295 c = gv*((
order_[p]-one)*p2*pf1*pf1 + p1*pf2);
306 std::vector<Real> &hvstat,
308 const std::vector<Real> &vstat,
310 const std::vector<Real> &xstat,
312 const Real
zero(0), one(1), two(2);
320 if ( pval_sum[p] >
zero ) {
323 c =
coeff_[p]*(pgv_sum[p]/std::pow(pval_sum[p],two-one/
order_[p]));
Provides the interface to evaluate objective functions.
MeanDeviationFromTarget(const Real target, const Real order, const Real coeff, const Ptr< PositiveFunction< Real > > &pf)
Constructor.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
void computeHessVec(Vector< Real > &hv, Objective< Real > &obj, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
typename PV< Real >::size_type size_type
std::vector< Ptr< Vector< Real > > > pg_
Real computeValue(Objective< Real > &obj, const Vector< Real > &x, Real &tol)
std::vector< Real > target_
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Ptr< Vector< Real > > hv_
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.
Ptr< PositiveFunction< Real > > positiveFunction_
void initialize(const Vector< Real > &x)
Initialize temporary variables.
MeanDeviationFromTarget(ROL::ParameterList &parlist)
Constructor.
Ptr< Vector< Real > > dualVector_
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 an interface for the mean plus a sum of arbitrary order deviations from targets.
Defines the linear algebra or vector space interface.
void updateValue(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal storage for value computation.
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.
void sumAll(Real *input, Real *output, int dim) const
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
std::vector< Real > coeff_
std::vector< Real >::size_type uint
std::vector< Ptr< Vector< Real > > > pg0_
std::vector< Real > pval_
void computeGradient(Vector< Real > &g, Objective< Real > &obj, const Vector< Real > &x, Real &tol)
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 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.
std::vector< Real > order_
std::vector< Ptr< Vector< Real > > > phv_
Real getValue(const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure value.
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.
MeanDeviationFromTarget(const std::vector< Real > &target, const std::vector< Real > &order, const std::vector< Real > &coeff, const Ptr< PositiveFunction< Real > > &pf)
Constructor.