44 #ifndef ROL_MEANVARIANCEFROMTARGET_HPP
45 #define ROL_MEANVARIANCEFROMTARGET_HPP
52 #include "ROL_ParameterList.hpp"
102 ROL_TEST_FOR_EXCEPTION((oSize!=cSize),std::invalid_argument,
103 ">>> ERROR (ROL::MeanVarianceFromTarget): Order and coefficient arrays have different sizes!");
104 Real
zero(0), two(2);
105 for (
int i = 0; i < oSize; i++) {
106 ROL_TEST_FOR_EXCEPTION((
order_[i] < two), std::invalid_argument,
107 ">>> ERROR (ROL::MeanVarianceFromTarget): Element of order array out of range!");
108 ROL_TEST_FOR_EXCEPTION((
coeff_[i] <
zero), std::invalid_argument,
109 ">>> ERROR (ROL::MeanVarianceFromTarget): Element of coefficient array out of range!");
112 ">>> ERROR (ROL::MeanVarianceFromTarget): PositiveFunction pointer is null!");
147 const std::vector<Real> &order,
148 const std::vector<Real> &coeff,
152 for (
uint i = 0; i < target.size(); i++ ) {
155 for (
uint i = 0; i < order.size(); i++ ) {
156 order_.push_back(order[i]);
158 for (
uint i = 0; i < coeff.size(); i++ ) {
159 coeff_.push_back(coeff[i]);
179 ROL::ParameterList &list
180 = parlist.sublist(
"SOL").sublist(
"Risk Measure").sublist(
"Mean Plus Variance From Target");
182 target_ = ROL::getArrayFromStringParameter<double>(list,
"Targets");
183 order_ = ROL::getArrayFromStringParameter<double>(list,
"Orders");
184 coeff_ = ROL::getArrayFromStringParameter<double>(list,
"Coefficients");
187 std::string type = list.get<std::string>(
"Deviation Type");
188 if ( type ==
"Upper" ) {
191 else if ( type ==
"Absolute" ) {
195 ROL_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
196 ">>> (ROL::MeanDeviation): Deviation type is not recoginized!");
205 const std::vector<Real> &xstat,
207 Real diff(0), pf0(0);
219 const std::vector<Real> &xstat,
221 Real diff(0), pf0(0), pf1(0), c(1), one(1);
235 const std::vector<Real> &vstat,
237 const std::vector<Real> &xstat,
239 Real diff(0), pf0(0), pf1(0), pf2(0), p1(0), p2(0), ch(1), cg(0), one(1), two(2);
248 p1 = std::pow(pf0,
order_[p]-one);
249 p2 = std::pow(pf0,
order_[p]-two);
Provides the interface to evaluate objective functions.
MeanVarianceFromTarget(const std::vector< Real > &target, const std::vector< Real > &order, const std::vector< Real > &coeff, const ROL::Ptr< PositiveFunction< Real > > &pf)
Constructor.
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
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.
Real computeValue(Objective< Real > &obj, const Vector< Real > &x, Real &tol)
Ptr< Vector< Real > > hv_
Ptr< Vector< Real > > dualVector_
Defines the linear algebra or vector space interface.
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
void updateValue(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal storage for value computation.
std::vector< Real > target_
void checkInputs(void) const
MeanVarianceFromTarget(ROL::ParameterList &parlist)
Constructor.
void updateGradient(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal risk measure storage for gradient computation.
std::vector< Real > order_
MeanVarianceFromTarget(const Real target, const Real order, const Real coeff, const ROL::Ptr< PositiveFunction< Real > > &pf)
Constructor.
void computeGradient(Vector< Real > &g, Objective< Real > &obj, const Vector< Real > &x, Real &tol)
Provides an interface for the mean plus a sum of arbitrary order variances from targets.
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...
std::vector< Real > coeff_
std::vector< Real >::size_type uint
ROL::Ptr< PositiveFunction< Real > > positiveFunction_