ROL
ROL_MeanDeviationFromTarget.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Rapid Optimization Library (ROL) Package
4 //
5 // Copyright 2014 NTESS and the ROL contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef ROL_MEANDEVIATIONFROMTARGET_HPP
11 #define ROL_MEANDEVIATIONFROMTARGET_HPP
12 
14 #include "ROL_ParameterList.hpp"
15 #include "ROL_PositiveFunction.hpp"
16 #include "ROL_PlusFunction.hpp"
17 #include "ROL_AbsoluteValue.hpp"
18 
39 namespace ROL {
40 
41 template<class Real>
44 private:
45  Ptr<PositiveFunction<Real> > positiveFunction_;
46  std::vector<Real> target_;
47  std::vector<Real> order_;
48  std::vector<Real> coeff_;
50 
51  std::vector<Real> pval_;
52  std::vector<Real> pgv_;
53 
54  std::vector<Ptr<Vector<Real> > > pg0_;
55  std::vector<Ptr<Vector<Real> > > pg_;
56  std::vector<Ptr<Vector<Real> > > phv_;
57 
59 
65 
68 
73 
74  void initializeMDT(void) {
75  // Initialize additional storage
76  pg_.clear(); pg0_.clear(); phv_.clear(); pval_.clear(); pgv_.clear();
77  pg_.resize(NumMoments_);
78  pg0_.resize(NumMoments_);
79  phv_.resize(NumMoments_);
80  pval_.resize(NumMoments_);
81  pgv_.resize(NumMoments_);
82  }
83 
84  void checkInputs(void) {
85  int oSize = order_.size(), cSize = coeff_.size(), tSize = target_.size();
86  ROL_TEST_FOR_EXCEPTION((oSize!=cSize),std::invalid_argument,
87  ">>> ERROR (ROL::MeanDeviationFromTarget): Order and coefficient arrays have different sizes!");
88  ROL_TEST_FOR_EXCEPTION((oSize!=tSize),std::invalid_argument,
89  ">>> ERROR (ROL::MeanDeviationFromTarget): Order and target arrays have different sizes!");
90  Real zero(0), two(2);
91  for (int i = 0; i < oSize; i++) {
92  ROL_TEST_FOR_EXCEPTION((order_[i] < two), std::invalid_argument,
93  ">>> ERROR (ROL::MeanDeviationFromTarget): Element of order array out of range!");
94  ROL_TEST_FOR_EXCEPTION((coeff_[i] < zero), std::invalid_argument,
95  ">>> ERROR (ROL::MeanDeviationFromTarget): Element of coefficient array out of range!");
96  }
97  ROL_TEST_FOR_EXCEPTION(positiveFunction_ == nullPtr, std::invalid_argument,
98  ">>> ERROR (ROL::MeanDeviationFromTarget): PositiveFunction pointer is null!");
99  initializeMDT();
100  }
101 
102 public:
113  MeanDeviationFromTarget( const Real target, const Real order, const Real coeff,
114  const Ptr<PositiveFunction<Real> > &pf )
115  : RandVarFunctional<Real>(), positiveFunction_(pf), firstResetMDT_(true) {
116  order_.clear(); order_.push_back(order);
117  coeff_.clear(); coeff_.push_back(coeff);
118  target_.clear(); target_.push_back(target);
119  NumMoments_ = order_.size();
120  checkInputs();
121  }
122 
133  MeanDeviationFromTarget( const std::vector<Real> &target,
134  const std::vector<Real> &order,
135  const std::vector<Real> &coeff,
136  const Ptr<PositiveFunction<Real> > &pf )
137  : RandVarFunctional<Real>(), positiveFunction_(pf), firstResetMDT_(true) {
138  target_.clear(); order_.clear(); coeff_.clear();
139  for ( uint i = 0; i < target.size(); i++ ) {
140  target_.push_back(target[i]);
141  }
142  for ( uint i = 0; i < order.size(); i++ ) {
143  order_.push_back(order[i]);
144  }
145  for ( uint i = 0; i < coeff.size(); i++ ) {
146  coeff_.push_back(coeff[i]);
147  }
148  NumMoments_ = order_.size();
149  checkInputs();
150  }
151 
164  MeanDeviationFromTarget( ROL::ParameterList &parlist )
165  : RandVarFunctional<Real>(), firstResetMDT_(true) {
166  ROL::ParameterList &list
167  = parlist.sublist("SOL").sublist("Risk Measure").sublist("Mean Plus Deviation From Target");
168 
169  // Get data from parameter list
170  target_ = ROL::getArrayFromStringParameter<double>(list,"Targets");
171 
172  order_ = ROL::getArrayFromStringParameter<double>(list,"Orders");
173 
174  coeff_ = ROL::getArrayFromStringParameter<double>(list,"Coefficients");
175 
176  // Build (approximate) positive function
177  std::string type = list.get<std::string>("Deviation Type");
178  if ( type == "Upper" ) {
179  positiveFunction_ = makePtr<PlusFunction<Real>>(list);
180  }
181  else if ( type == "Absolute" ) {
182  positiveFunction_ = makePtr<AbsoluteValue<Real>>(list);
183  }
184  else {
185  ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
186  ">>> (ROL::MeanDeviation): Deviation type is not recoginized!");
187  }
188  // Check inputs
189  NumMoments_ = order_.size();
190  checkInputs();
191  }
192 
193  void initialize(const Vector<Real> &x) {
195  if (firstResetMDT_) {
196  for ( uint p = 0; p < NumMoments_; p++ ) {
197  pg0_[p] = x.dual().clone();
198  pg_[p] = x.dual().clone();
199  phv_[p] = x.dual().clone();
200  }
201  firstResetMDT_ = false;
202  }
203  Real zero(0);
204  for ( uint p = 0; p < NumMoments_; p++ ) {
205  pg0_[p]->zero(); pg_[p]->zero(); phv_[p]->zero();
206  pval_[p] = zero; pgv_[p] = zero;
207  }
208  }
209 
211  const Vector<Real> &x,
212  const std::vector<Real> &xstat,
213  Real &tol) {
214  Real diff(0), pf0(0);
215  Real val = computeValue(obj,x,tol);
216  val_ += weight_ * val;
217  for ( uint p = 0; p < NumMoments_; p++ ) {
218  diff = val-target_[p];
219  pf0 = positiveFunction_->evaluate(diff,0);
220  pval_[p] += weight_ * std::pow(pf0,order_[p]);
221  }
222  }
223 
224  Real getValue(const Vector<Real> &x,
225  const std::vector<Real> &xstat,
226  SampleGenerator<Real> &sampler) {
227  const Real one(1);
228  Real dev(0);
229  sampler.sumAll(&val_,&dev,1);
230  std::vector<Real> pval_sum(NumMoments_);
231  sampler.sumAll(&pval_[0],&pval_sum[0],NumMoments_);
232  for ( uint p = 0; p < NumMoments_; p++ ) {
233  dev += coeff_[p] * std::pow(pval_sum[p],one/order_[p]);
234  }
235  return dev;
236  }
237 
239  const Vector<Real> &x,
240  const std::vector<Real> &xstat,
241  Real &tol) {
242  Real diff(0), pf0(0), pf1(0), c(0), one(1);
243  Real val = computeValue(obj,x,tol);
244  computeGradient(*dualVector_,obj,x,tol);
245  for ( uint p = 0; p < NumMoments_; p++ ) {
246  diff = val-target_[p];
247  pf0 = positiveFunction_->evaluate(diff,0);
248  pf1 = positiveFunction_->evaluate(diff,1);
249  c = std::pow(pf0,order_[p]-one) * pf1;
250  (pg_[p])->axpy(weight_ * c,*dualVector_);
251  pval_[p] += weight_ * std::pow(pf0,order_[p]);
252  }
253  g_->axpy(weight_,*dualVector_);
254  }
255 
257  std::vector<Real> &gstat,
258  const Vector<Real> &x,
259  const std::vector<Real> &xstat,
260  SampleGenerator<Real> &sampler) {
261  const Real zero(0), one(1);
262  sampler.sumAll(*g_,g);
263  std::vector<Real> pval_sum(NumMoments_);
264  sampler.sumAll(&pval_[0],&pval_sum[0],NumMoments_);
265  for ( uint p = 0; p < NumMoments_; p++ ) {
266  if ( pval_sum[p] > zero ) {
267  dualVector_->zero();
268  sampler.sumAll(*(pg_[p]),*dualVector_);
269  g.axpy(coeff_[p]/std::pow(pval_sum[p],one-one/order_[p]),*dualVector_);
270  }
271  }
272  }
273 
275  const Vector<Real> &v,
276  const std::vector<Real> &vstat,
277  const Vector<Real> &x,
278  const std::vector<Real> &xstat,
279  Real &tol) {
280  const Real one(1), two(2);
281  Real diff(0), pf0(0), pf1(0), pf2(0), p0(0), p1(0), p2(0), c(0);
282  Real val = computeValue(obj,x,tol);
283  Real gv = computeGradVec(*g_,obj,v,x,tol);
284  computeHessVec(*dualVector_,obj,v,x,tol);
285  for ( uint p = 0; p < NumMoments_; p++ ) {
286  diff = val - target_[p];
287  pf0 = positiveFunction_->evaluate(diff,0);
288  pf1 = positiveFunction_->evaluate(diff,1);
289  pf2 = positiveFunction_->evaluate(diff,2);
290  p0 = std::pow(pf0,order_[p]);
291  p1 = std::pow(pf0,order_[p]-one);
292  p2 = std::pow(pf0,order_[p]-two);
293  c = -(order_[p]-one)*p1*pf1;
294  pg0_[p]->axpy(weight_*c,*g_);
295  c = gv*((order_[p]-one)*p2*pf1*pf1 + p1*pf2);
296  pg_[p]->axpy(weight_*c,*g_);
297  c = p1*pf1;
298  phv_[p]->axpy(weight_*c,*dualVector_);
299  pval_[p] += weight_*p0;
300  pgv_[p] += weight_*p1*pf1*gv;
301  }
302  hv_->axpy(weight_,*dualVector_);
303  }
304 
306  std::vector<Real> &hvstat,
307  const Vector<Real> &v,
308  const std::vector<Real> &vstat,
309  const Vector<Real> &x,
310  const std::vector<Real> &xstat,
311  SampleGenerator<Real> &sampler) {
312  const Real zero(0), one(1), two(2);
313  sampler.sumAll(*hv_,hv);
314  std::vector<Real> pval_sum(NumMoments_);
315  sampler.sumAll(&(pval_)[0],&pval_sum[0],NumMoments_);
316  std::vector<Real> pgv_sum(NumMoments_);
317  sampler.sumAll(&(pgv_)[0],&pgv_sum[0],NumMoments_);
318  Real c(0);
319  for ( uint p = 0; p < NumMoments_; p++ ) {
320  if ( pval_sum[p] > zero ) {
321  dualVector_->zero();
322  sampler.sumAll(*(pg0_[p]),*dualVector_);
323  c = coeff_[p]*(pgv_sum[p]/std::pow(pval_sum[p],two-one/order_[p]));
324  hv.axpy(c,*dualVector_);
325 
326  dualVector_->zero();
327  sampler.sumAll(*(pg_[p]),*dualVector_);
328  c = coeff_[p]/std::pow(pval_sum[p],one-one/order_[p]);
329  hv.axpy(c,*dualVector_);
330 
331  dualVector_->zero();
332  sampler.sumAll(*(phv_[p]),*dualVector_);
333  hv.axpy(c,*dualVector_);
334  }
335  }
336  }
337 };
338 
339 }
340 
341 #endif
Provides the interface to evaluate objective functions.
MeanDeviationFromTarget(const Real target, const Real order, const Real coeff, const Ptr< PositiveFunction< Real > > &pf)
Constructor.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: ROL_Vector.hpp:192
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
Ptr< Vector< Real > > g_
std::vector< Ptr< Vector< Real > > > pg_
Real computeValue(Objective< Real > &obj, const Vector< Real > &x, Real &tol)
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:119
Ptr< Vector< Real > > hv_
void getHessVec(Vector< Real > &hv, std::vector< Real > &hvstat, const Vector< Real > &v, const std::vector< Real > &vstat, const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure Hessian-times-a-vector.
Ptr< PositiveFunction< Real > > positiveFunction_
void initialize(const Vector< Real > &x)
Initialize temporary variables.
MeanDeviationFromTarget(ROL::ParameterList &parlist)
Constructor.
Ptr< Vector< Real > > dualVector_
void updateGradient(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal risk measure storage for gradient computation.
Provides an interface for the mean plus a sum of arbitrary order deviations from targets.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:46
void updateValue(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal storage for value computation.
void getGradient(Vector< Real > &g, std::vector< Real > &gstat, const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure (sub)gradient.
void sumAll(Real *input, Real *output, int dim) const
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
std::vector< Real >::size_type uint
std::vector< Ptr< Vector< Real > > > pg0_
void computeGradient(Vector< Real > &g, Objective< Real > &obj, const Vector< Real > &x, Real &tol)
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...
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.
std::vector< Ptr< Vector< Real > > > phv_
Real getValue(const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure value.
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.
MeanDeviationFromTarget(const std::vector< Real > &target, const std::vector< Real > &order, const std::vector< Real > &coeff, const Ptr< PositiveFunction< Real > > &pf)
Constructor.