44 #ifndef ROL_RISKNEUTRALOBJECTIVE_HPP
45 #define ROL_RISKNEUTRALOBJECTIVE_HPP
82 const std::vector<Real> ¶m, Real &tol) {
90 value_storage_.insert(std::pair<std::vector<Real>,Real>(param,val));
96 const std::vector<Real> ¶m, Real &tol) {
104 Ptr<Vector<Real>> tmp = g.
clone();
112 const std::vector<Real> ¶m, Real &tol) {
123 const bool storage =
true )
134 const bool storage =
true )
144 const bool storage =
true )
157 value_ =
static_cast<Real
>(0);
173 Real myval(0), ptval(0), val(0), one(1), two(2), error(two*tol + one);
174 std::vector<Real> ptvals;
175 while ( error > tol ) {
180 ptvals.push_back(ptval);
195 std::vector<Ptr<Vector<Real>>> ptgs;
196 Real one(1), two(2), error(two*tol + one);
197 while ( error > tol ) {
Provides the interface to evaluate objective functions.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
RiskNeutralObjective(const Ptr< Objective< Real >> &pObj, const Ptr< SampleGenerator< Real >> &vsampler, const Ptr< SampleGenerator< Real >> &gsampler, const Ptr< SampleGenerator< Real >> &hsampler, const bool storage=true)
Ptr< Vector< Real > > pointDual_
Ptr< Vector< Real > > gradient_
void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply preconditioner to vector.
RiskNeutralObjective(const Ptr< Objective< Real >> &pObj, const Ptr< SampleGenerator< Real >> &sampler, const bool storage=true)
void initialize(const Vector< Real > &x)
void getHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, const std::vector< Real > ¶m, Real &tol)
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
virtual void zero()
Set to zero vector.
Ptr< SampleGenerator< Real > > ValueSampler_
Defines the linear algebra or vector space interface.
RiskNeutralObjective(const Ptr< Objective< Real >> &pObj, const Ptr< SampleGenerator< Real >> &vsampler, const Ptr< SampleGenerator< Real >> &gsampler, const bool storage=true)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
void getValue(Real &val, const Vector< Real > &x, const std::vector< Real > ¶m, Real &tol)
Ptr< Objective< Real > > ParametrizedObjective_
Ptr< Vector< Real > > sumDual_
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
Ptr< SampleGenerator< Real > > GradientSampler_
std::map< std::vector< Real >, Ptr< Vector< Real > > > gradient_storage_
std::map< std::vector< Real >, Real > value_storage_
void getGradient(Vector< Real > &g, const Vector< Real > &x, const std::vector< Real > ¶m, Real &tol)
virtual void set(const Vector &x)
Set where .
Ptr< SampleGenerator< Real > > HessianSampler_