44 #ifndef PH_BPOEOBJECTIVE_H
45 #define PH_BPOEOBJECTIVE_H
60 const Ptr<Objective<Real>>
obj_;
101 Ptr<const std::vector<Real>> xstat = xrv.
getStatistic();
102 if (xstat == nullPtr) {
103 xstat = makePtr<const std::vector<Real>>(0);
111 if (xstat == nullPtr) {
112 xstat = makePtr<std::vector<Real>>(0);
118 Real
pplus(
const Real x,
const int deriv = 0)
const {
119 const Real
zero(0), one(1), two(2), three(3);
129 : std::pow(x,order_-one)));
135 : std::pow(x,
order_-two))));
144 return pplus(arg,deriv);
150 ParameterList &parlist)
155 ParameterList &list = parlist.sublist(
"SOL").sublist(
"Probability").sublist(
"bPOE");
157 order_ = list.get<Real>(
"Moment Order");
162 obj_->update(*xvec,flag,iter);
171 Real xt = (*xstat)[0];
178 Ptr<std::vector<Real>> gstat =
getStat(g);
186 gvec->scale(prob*xt);
187 (*gstat)[0] = prob*diff;
192 Ptr<std::vector<Real>> hstat =
getStat(hv);
202 Real gv = vvec->dot(
g_->dual());
203 obj_->hessVec(*hvec,*vvec,*xvec,tol);
204 hvec->scale(prob1*xt);
205 hvec->axpy(prob2*xt*(vt*diff+xt*gv)+vt*prob1,*
g_);
206 (*hstat)[0] = prob2*std::pow(diff,2)*vt+(prob2*diff*xt+prob1)*gv;
210 obj_->setParameter(param);
bool isGradientInitialized_
Provides the interface to evaluate objective functions.
ROL::Ptr< std::vector< Real > > getStatistic(const int comp=0, const int index=0)
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
void getValue(const Vector< Real > &x, Real &tol)
ROL::Ptr< const Vector< Real > > getVector(void) const
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
Real pplus(const Real x, const int deriv=0) const
void setParameter(const std::vector< Real > ¶m)
Defines the linear algebra or vector space interface.
Provides the interface for the progressive hedging probability objective.
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Ptr< Vector< Real > > getVector(Vector< Real > &x) const
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
Ptr< std::vector< Real > > getStat(Vector< Real > &x) const
PH_bPOEObjective(const Ptr< Objective< Real >> &obj, ParameterList &parlist)
Real value(const Vector< Real > &x, Real &tol)
Compute value.
void getGradient(const Vector< Real > &x, Real &tol)
Real bPOEobjective(const Real t, const Real x, const int deriv=0) const
virtual void setParameter(const std::vector< Real > ¶m)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Ptr< const Vector< Real > > getConstVector(const Vector< Real > &x) const
const Ptr< Objective< Real > > obj_
Ptr< const std::vector< Real > > getConstStat(const Vector< Real > &x) const