44 #ifndef ROL_HMCROBJECTIVE_HPP 45 #define ROL_HMCROBJECTIVE_HPP 47 #include "Teuchos_RCP.hpp" 85 pointGrad_ = x.
dual().clone();
86 pointHess_ = x.
dual().clone();
87 gradient0_ = x.
dual().clone();
88 sumGrad0_ = x.
dual().clone();
89 gradient1_ = x.
dual().clone();
90 sumGrad1_ = x.
dual().clone();
91 gradient2_ = x.
dual().clone();
92 sumGrad2_ = x.
dual().clone();
93 hessvec_ = x.
dual().clone();
94 sumHess_ = x.
dual().clone();
102 if ( !initialized_ ) {
108 const std::vector<Real> ¶m, Real &tol) {
109 if ( storage_ && value_storage_.count(param) ) {
110 val = value_storage_[param];
113 ParametrizedObjective_->setParameter(param);
114 val = ParametrizedObjective_->value(x,tol);
116 value_storage_.insert(std::pair<std::vector<Real>,Real>(param,val));
122 const std::vector<Real> ¶m, Real &tol) {
123 if ( storage_ && gradient_storage_.count(param) ) {
124 g.
set(*(gradient_storage_[param]));
127 ParametrizedObjective_->setParameter(param);
128 ParametrizedObjective_->gradient(g,x,tol);
130 Teuchos::RCP<Vector<Real> > tmp = g.
clone();
131 gradient_storage_.insert(std::pair<std::vector<Real>,Teuchos::RCP<
Vector<Real> > >(param,tmp));
132 gradient_storage_[param]->set(g);
138 const std::vector<Real> ¶m, Real &tol) {
139 ParametrizedObjective_->setParameter(param);
140 ParametrizedObjective_->hessVec(hv,v,x,tol);
148 Real order, Real prob,
152 bool storage =
true )
153 : ParametrizedObjective_(pObj),
154 ValueSampler_(vsampler), GradientSampler_(gsampler), HessianSampler_(hsampler),
155 initialized_(false), storage_(storage) {
156 order_ = ((order < 1.0) ? 1.0 : order);
157 prob_ = ((prob > 1.0) ? 1.0 : ((prob < 0.0) ? 0.0 : prob));
158 value_storage_.clear();
159 gradient_storage_.clear();
163 Real order, Real prob,
166 bool storage =
true )
167 : ParametrizedObjective_(pObj),
168 ValueSampler_(vsampler), GradientSampler_(gsampler), HessianSampler_(gsampler),
169 initialized_(false), storage_(storage) {
170 order_ = ((order < 1.0) ? 1.0 : order);
171 prob_ = ((prob > 1.0) ? 1.0 : ((prob < 0.0) ? 0.0 : prob));
172 value_storage_.clear();
173 gradient_storage_.clear();
177 Real order, Real prob,
179 bool storage =
true )
180 : ParametrizedObjective_(pObj), order_(order), prob_(prob),
181 ValueSampler_(sampler), GradientSampler_(sampler), HessianSampler_(sampler),
182 initialized_(false), storage_(storage) {
183 order_ = ((order < 1.0) ? 1.0 : order);
184 prob_ = ((prob > 1.0) ? 1.0 : ((prob < 0.0) ? 0.0 : prob));
185 value_storage_.clear();
186 gradient_storage_.clear();
190 Teuchos::RCP<Vector<Real> > xvec; Real xvar = 0.0;
192 ParametrizedObjective_->update(*xvec,flag,iter);
193 ValueSampler_->update(*xvec);
195 value_storage_.clear();
198 GradientSampler_->update(*xvec);
199 HessianSampler_->update(*xvec);
201 gradient_storage_.clear();
207 Teuchos::RCP<Vector<Real> > xvec; Real xvar = 0.0;
210 std::vector<Real> point;
211 Real weight = 0.0, myval = 0.0, pval = 0.0, val = 0.0;
212 int start = ValueSampler_->start(), end = ValueSampler_->numMySamples();
213 for (
int i = start; i < end; i++ ) {
214 weight = ValueSampler_->getMyWeight(i);
215 point = ValueSampler_->getMyPoint(i);
220 myval += weight*((order_ == 1.0) ? pval-xvar
221 : std::pow(pval-xvar,order_));
225 ValueSampler_->sumAll(&myval,&val,1);
227 if (std::abs(val) < ROL_EPSILON<Real>()) {
230 return xvar + ((order_ == 1.0) ? val
231 : ((order_ == 2.0) ? std::sqrt(val)
232 : std::pow(val,1.0/order_)))/(1.0-
prob_);
236 Teuchos::RCP<Vector<Real> > xvec; Real xvar = 0.0;
240 g.
zero(); sumGrad0_->zero();
241 std::vector<Real> point, val(2,0.0), myval(2,0.0);
242 Real weight = 0.0, pval = 0.0, pvalp0 = 0.0, pvalp1 = 0.0;
243 int start = GradientSampler_->start(), end = GradientSampler_->numMySamples();
244 for (
int i = start; i < end; i++ ) {
245 weight = GradientSampler_->getMyWeight(i);
246 point = GradientSampler_->getMyPoint(i);
251 pvalp0 = ((order_ == 1.0) ? pval-xvar
252 : std::pow(pval-xvar,order_));
253 pvalp1 = ((order_ == 1.0) ? 1.0
254 : ((order_ == 2.0) ? pval-xvar
255 : std::pow(pval-xvar,order_-1.0)));
257 myval[0] += weight*pvalp0;
258 myval[1] += weight*pvalp1;
262 sumGrad0_->axpy(weight*pvalp1,*pointGrad_);
265 Real gvar = 1.0; gradient0_->zero();
267 GradientSampler_->sumAll(&myval[0],&val[0],2);
268 if (std::abs(val[0]) >= ROL_EPSILON<Real>()) {
269 GradientSampler_->sumAll(*sumGrad0_,*gradient0_);
271 Real norm = ((order_ == 1.0) ? 1.0
272 : ((order_ == 2.0) ? std::sqrt(val[0])
273 : std::pow(val[0],(order_-1.0)/order_)));
274 gvar -= val[1]/((1.0-
prob_)*norm);
275 gradient0_->scale(1.0/((1.0-prob_)*norm));
284 Teuchos::RCP<Vector<Real> > xvec; Real xvar = 0.0;
286 Teuchos::RCP<Vector<Real> > vvec; Real vvar = 0.0;
291 sumGrad0_->zero(); sumGrad1_->zero(); sumGrad2_->zero(); sumHess_->zero();
292 gradient0_->zero(); gradient1_->zero(); gradient2_->zero();
294 std::vector<Real> point;
295 Real pval = 0.0, pvalp0 = 0.0, pvalp1 = 0.0, pvalp2 = 0.0, gv = 0.0;
296 std::vector<Real> val(5,0.0), myval(5,0.0);
297 int start = HessianSampler_->start(), end = HessianSampler_->numMySamples();
298 for (
int i = start; i < end; i++ ) {
300 weight = HessianSampler_->getMyWeight(i);
301 point = HessianSampler_->getMyPoint(i);
306 pvalp0 = ((order_ == 1.0) ? pval-xvar
307 : std::pow(pval-xvar,order_));
308 pvalp1 = ((order_ == 1.0) ? 1.0
309 : ((order_ == 2.0) ? pval-xvar
310 : std::pow(pval-xvar,order_-1.0)));
311 pvalp2 = ((order_ == 1.0) ? 0.0
312 : ((order_ == 2.0) ? 1.0
313 : ((order_ == 3.0) ? pval-xvar
314 : std::pow(pval-xvar,order_-2.0))));
316 myval[0] += weight*pvalp0;
317 myval[1] += weight*pvalp1;
318 myval[2] += weight*pvalp2;
321 gv = pointGrad_->dot(vvec->dual());
323 myval[3] += weight*pvalp1*gv;
324 myval[4] += weight*pvalp2*gv;
325 sumGrad0_->axpy(weight*pvalp1,*pointGrad_);
326 sumGrad1_->axpy(weight*pvalp2,*pointGrad_);
327 sumGrad2_->axpy(weight*pvalp2*gv,*pointGrad_);
329 getHessVec(*pointHess_,*vvec,*xvec,point,tol);
331 sumHess_->axpy(weight*pvalp1,*pointHess_);
334 Real hvar = 0.0; hessvec_->zero();
335 HessianSampler_->sumAll(&myval[0],&val[0],5);
336 if (std::abs(val[0]) >= ROL_EPSILON<Real>()) {
338 HessianSampler_->sumAll(*sumGrad0_,*gradient0_);
339 HessianSampler_->sumAll(*sumGrad1_,*gradient1_);
340 HessianSampler_->sumAll(*sumGrad2_,*gradient2_);
341 HessianSampler_->sumAll(*sumHess_,*hessvec_);
343 Real norm0 = (1.0-
prob_)*((order_ == 1.0) ? 1.0
344 : ((order_ == 2.0) ? std::sqrt(val[0])
345 : std::pow(val[0],(order_-1.0)/order_)));
346 Real norm1 = (1.0-
prob_)*((order_ == 1.0) ? 1.0
347 : std::pow(val[0],(2.0*order_-1.0)/order_));
348 hvar = (order_-1.0)*((val[2]/norm0 - val[1]*val[1]/norm1)*vvar
349 -(val[4]/norm0 - val[3]*val[1]/norm1));
350 hessvec_->scale(1.0/norm0);
351 hessvec_->axpy(-(order_-1.0)*vvar/norm0,*gradient1_);
352 hessvec_->axpy((order_-1.0)*(vvar*val[1]-val[3])/norm1,*gradient0_);
353 hessvec_->axpy((order_-1.0)/norm0,*gradient2_);
void getValue(Real &val, const Vector< Real > &x, const std::vector< Real > ¶m, Real &tol)
Provides the interface to evaluate objective functions.
Teuchos::RCP< SampleGenerator< Real > > GradientSampler_
Teuchos::RCP< SampleGenerator< Real > > ValueSampler_
void getHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, const std::vector< Real > ¶m, Real &tol)
Teuchos::RCP< Vector< Real > > sumGrad1_
Teuchos::RCP< Vector< Real > > gradient1_
Teuchos::RCP< Vector< Real > > hessvec_
Teuchos::RCP< Vector< Real > > sumGrad2_
Teuchos::RCP< Vector< Real > > gradient0_
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply preconditioner to vector.
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
Teuchos::RCP< Vector< Real > > pointGrad_
std::map< std::vector< Real >, Real > value_storage_
Teuchos::RCP< SampleGenerator< Real > > HessianSampler_
std::map< std::vector< Real >, Teuchos::RCP< Vector< Real > > > gradient_storage_
void setVector(const Vector< Real > &vec)
void initialize(const Vector< Real > &x)
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
void setStatistic(const Real stat)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
HMCRObjective(Teuchos::RCP< Objective< Real > > &pObj, Real order, Real prob, Teuchos::RCP< SampleGenerator< Real > > &sampler, bool storage=true)
Teuchos::RCP< Vector< Real > > sumHess_
HMCRObjective(Teuchos::RCP< Objective< Real > > &pObj, Real order, Real prob, Teuchos::RCP< SampleGenerator< Real > > &vsampler, Teuchos::RCP< SampleGenerator< Real > > &gsampler, Teuchos::RCP< SampleGenerator< Real > > &hsampler, bool storage=true)
Teuchos::RCP< Objective< Real > > ParametrizedObjective_
Teuchos::RCP< Vector< Real > > pointHess_
virtual void set(const Vector &x)
Set where .
HMCRObjective(Teuchos::RCP< Objective< Real > > &pObj, Real order, Real prob, Teuchos::RCP< SampleGenerator< Real > > &vsampler, Teuchos::RCP< SampleGenerator< Real > > &gsampler, bool storage=true)
Real value(const Vector< Real > &x, Real &tol)
Compute value.
Teuchos::RCP< Vector< Real > > gradient2_
void getGradient(Vector< Real > &g, const Vector< Real > &x, const std::vector< Real > ¶m, Real &tol)
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
void unwrap_const_CVaR_vector(Teuchos::RCP< Vector< Real > > &xvec, Real &xvar, const Vector< Real > &x)
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
Teuchos::RCP< Vector< Real > > sumGrad0_