44 #ifndef ROL_RISKNEUTRALOBJECTIVE_HPP 45 #define ROL_RISKNEUTRALOBJECTIVE_HPP 47 #include "Teuchos_RefCountPtr.hpp" 74 const std::vector<Real> ¶m, Real &tol) {
75 if ( storage_ && value_storage_.count(param) ) {
76 val = value_storage_[param];
79 ParametrizedObjective_->setParameter(param);
80 val = ParametrizedObjective_->value(x,tol);
82 value_storage_.insert(std::pair<std::vector<Real>,Real>(param,val));
88 const std::vector<Real> ¶m, Real &tol) {
89 if ( storage_ && gradient_storage_.count(param) ) {
90 g.
set(*(gradient_storage_[param]));
93 ParametrizedObjective_->setParameter(param);
94 ParametrizedObjective_->gradient(g,x,tol);
96 Teuchos::RCP<Vector<Real> > tmp = g.
clone();
97 gradient_storage_.insert(std::pair<std::vector<Real>,Teuchos::RCP<
Vector<Real> > >(param,tmp));
98 gradient_storage_[param]->set(g);
104 const std::vector<Real> ¶m, Real &tol) {
105 ParametrizedObjective_->setParameter(param);
106 ParametrizedObjective_->hessVec(hv,v,x,tol);
117 const bool storage =
true )
118 : ParametrizedObjective_(pObj),
119 ValueSampler_(vsampler), GradientSampler_(gsampler), HessianSampler_(hsampler),
120 firstUpdate_(true), storage_(storage) {
121 value_storage_.clear();
122 gradient_storage_.clear();
128 const bool storage =
true )
129 : ParametrizedObjective_(pObj),
130 ValueSampler_(vsampler), GradientSampler_(gsampler), HessianSampler_(gsampler),
131 firstUpdate_(true), storage_(storage) {
132 value_storage_.clear();
133 gradient_storage_.clear();
138 const bool storage =
true )
139 : ParametrizedObjective_(pObj),
140 ValueSampler_(sampler), GradientSampler_(sampler), HessianSampler_(sampler),
141 firstUpdate_(true), storage_(storage) {
142 value_storage_.clear();
143 gradient_storage_.clear();
147 if ( firstUpdate_ ) {
148 gradient_ = (x.
dual()).clone();
149 pointDual_ = (x.
dual()).clone();
150 sumDual_ = (x.
dual()).clone();
151 firstUpdate_ =
false;
153 ParametrizedObjective_->update(x,(flag && iter>=0),iter);
154 ValueSampler_->update(x);
155 value_ =
static_cast<Real
>(0);
157 value_storage_.clear();
159 if ( flag && iter>=0 ) {
160 GradientSampler_->update(x);
161 HessianSampler_->update(x);
164 gradient_storage_.clear();
170 Real myval(0), ptval(0), val(0), one(1), two(2), error(two*tol + one);
171 std::vector<Real> ptvals;
172 while ( error > tol ) {
173 ValueSampler_->refine();
174 for (
int i = ValueSampler_->start(); i < ValueSampler_->numMySamples(); ++i ) {
175 getValue(ptval,x,ValueSampler_->getMyPoint(i),tol);
176 myval += ValueSampler_->getMyWeight(i)*ptval;
177 ptvals.push_back(ptval);
179 error = ValueSampler_->computeError(ptvals);
182 ValueSampler_->sumAll(&myval,&val,1);
184 ValueSampler_->setSamples();
190 g.
zero(); pointDual_->zero(); sumDual_->zero();
191 std::vector<Teuchos::RCP<Vector<Real> > > ptgs;
192 Real one(1), two(2), error(two*tol + one);
193 while ( error > tol ) {
194 GradientSampler_->refine();
195 for (
int i = GradientSampler_->start(); i < GradientSampler_->numMySamples(); ++i ) {
196 getGradient(*pointDual_,x,GradientSampler_->getMyPoint(i),tol);
197 sumDual_->axpy(GradientSampler_->getMyWeight(i),*
pointDual_);
198 ptgs.push_back(pointDual_->clone());
199 (ptgs.back())->
set(*pointDual_);
201 error = GradientSampler_->computeError(ptgs,x);
207 GradientSampler_->sumAll(*sumDual_,g);
210 GradientSampler_->setSamples();
216 hv.
zero(); pointDual_->zero(); sumDual_->zero();
217 for (
int i = 0; i < HessianSampler_->numMySamples(); ++i ) {
218 getHessVec(*pointDual_,v,x,HessianSampler_->getMyPoint(i),tol);
219 sumDual_->axpy(HessianSampler_->getMyWeight(i),*
pointDual_);
221 HessianSampler_->sumAll(*sumDual_,hv);
Provides the interface to evaluate objective functions.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
RiskNeutralObjective(const Teuchos::RCP< Objective< Real > > &pObj, const Teuchos::RCP< SampleGenerator< Real > > &vsampler, const Teuchos::RCP< SampleGenerator< Real > > &gsampler, const Teuchos::RCP< SampleGenerator< Real > > &hsampler, const bool storage=true)
Teuchos::RCP< SampleGenerator< Real > > GradientSampler_
Teuchos::RCP< Vector< Real > > sumDual_
Teuchos::RCP< SampleGenerator< Real > > ValueSampler_
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
void getHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, const std::vector< Real > ¶m, Real &tol)
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
Teuchos::RCP< Vector< Real > > gradient_
Teuchos::RCP< Objective< Real > > ParametrizedObjective_
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual ~RiskNeutralObjective()
Teuchos::RCP< Vector< Real > > pointDual_
virtual 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)
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply preconditioner to vector.
std::map< std::vector< Real >, Real > value_storage_
std::map< std::vector< Real >, Teuchos::RCP< Vector< Real > > > gradient_storage_
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
void getGradient(Vector< Real > &g, const Vector< Real > &x, const std::vector< Real > ¶m, Real &tol)
RiskNeutralObjective(const Teuchos::RCP< Objective< Real > > &pObj, const Teuchos::RCP< SampleGenerator< Real > > &sampler, const bool storage=true)
RiskNeutralObjective(const Teuchos::RCP< Objective< Real > > &pObj, const Teuchos::RCP< SampleGenerator< Real > > &vsampler, const Teuchos::RCP< SampleGenerator< Real > > &gsampler, const bool storage=true)
virtual Real value(const Vector< Real > &x, Real &tol)
Compute value.
virtual void set(const Vector &x)
Set where .
Teuchos::RCP< SampleGenerator< Real > > HessianSampler_