44 #ifndef ROL_MEANDEVIATIONFROMTARGET_HPP 45 #define ROL_MEANDEVIATIONFROMTARGET_HPP 63 std::vector<Teuchos::RCP<Vector<Real> > >
pg0_;
64 std::vector<Teuchos::RCP<Vector<Real> > >
pg_;
65 std::vector<Teuchos::RCP<Vector<Real> > >
phv_;
73 :
RiskMeasure<Real>(), positiveFunction_(pf), firstReset_(true) {
75 target_.clear(); order_.clear(); coeff_.clear();
76 target_.push_back(target);
77 order_.push_back((order < 2.0) ? 2.0 : order);
78 coeff_.push_back((coeff < 0.0) ? 1.0 : coeff);
80 pg_.clear(); pg0_.clear(); phv_.clear(); pval_.clear(); pgv_.clear();
81 pg_.resize(order_.size());
82 pg0_.resize(order_.size());
83 phv_.resize(order_.size());
84 pval_.resize(order_.size(),0.0);
85 pgv_.resize(order_.size(),0.0);
90 :
RiskMeasure<Real>(), positiveFunction_(pf), firstReset_(true) {
92 target_.clear(); order_.clear(); coeff_.clear();
93 if ( order.size() != target.size() ) {
94 target.resize(order.size(),0.0);
96 if ( order.size() != coeff.size() ) {
97 coeff.resize(order.size(),1.0);
99 for (
unsigned i = 0; i < order.size(); i++ ) {
100 target_.push_back(target[i]);
101 order_.push_back((order[i] < 2.0) ? 2.0 : order[i]);
102 coeff_.push_back((coeff[i] < 0.0) ? 1.0 : coeff[i]);
105 pg_.clear(); pg0_.clear(); phv_.clear(); pval_.clear(); pgv_.clear();
106 pg_.resize(order_.size());
107 pg0_.resize(order_.size());
108 phv_.resize(order_.size());
109 pval_.resize(order_.size(),0.0);
110 pgv_.resize(order_.size(),0.0);
116 for (
unsigned p = 0; p < order_.size(); p++ ) {
117 pg0_[p] = (x.
dual()).clone();
118 pg_[p] = (x.
dual()).clone();
119 phv_[p] = (x.
dual()).clone();
123 for (
unsigned p = 0; p < order_.size(); p++ ) {
124 pg0_[p]->zero(); pg_[p]->zero(); phv_[p]->zero();
125 pval_[p] = 0.0; pgv_[p] = 0.0;
133 for (
unsigned p = 0; p < order_.size(); p++ ) {
134 pg0_[p] = (x.
dual()).clone();
135 pg_[p] = (x.
dual()).clone();
136 phv_[p] = (x.
dual()).clone();
140 for (
unsigned p = 0; p < order_.size(); p++ ) {
141 pg0_[p]->zero(); pg_[p]->zero(); phv_[p]->zero();
142 pval_[p] = 0.0; pgv_[p] = 0.0;
146 void update(
const Real val,
const Real weight) {
147 Real diff = 0.0, pf0 = 0.0;
149 for (
unsigned p = 0; p < order_.size(); p++ ) {
150 diff = val-target_[p];
151 pf0 = positiveFunction_->evaluate(diff,0);
152 pval_[p] += weight * std::pow(pf0,order_[p]);
157 Real diff = 0.0, pf0 = 0.0, pf1 = 0.0, c = 0.0;
158 for (
unsigned p = 0; p < order_.size(); p++ ) {
159 diff = val-target_[p];
160 pf0 = positiveFunction_->evaluate(diff,0);
161 pf1 = positiveFunction_->evaluate(diff,1);
162 c = std::pow(pf0,order_[p]-1.0) * pf1;
163 (pg_[p])->axpy(weight * c,g);
164 pval_[p] += weight * std::pow(pf0,order_[p]);
171 Real diff = 0.0, pf0 = 0.0, pf1 = 0.0, pf2 = 0.0, p0 = 0.0, p1 = 0.0, p2 = 0.0, c = 0.0;
172 for (
unsigned p = 0; p < order_.size(); p++ ) {
173 diff = val - target_[p];
174 pf0 = positiveFunction_->evaluate(diff,0);
175 pf1 = positiveFunction_->evaluate(diff,1);
176 pf2 = positiveFunction_->evaluate(diff,2);
177 p0 = std::pow(pf0,order_[p]);
178 p1 = std::pow(pf0,order_[p]-1.0);
179 p2 = std::pow(pf0,order_[p]-2.0);
180 c = -(order_[p]-1.0)*p1*pf1;
181 pg0_[p]->axpy(weight*c,g);
182 c = gv*((order_[p]-1.0)*p2*pf1*pf1 + p1*pf2);
183 pg_[p]->axpy(weight*c,g);
185 phv_[p]->axpy(weight*c,hv);
186 pval_[p] += weight*p0;
187 pgv_[p] += weight*p1*pf1*gv;
195 sampler.
sumAll(&val,&dev,1);
196 std::vector<Real> pval_sum(pval_.size());
197 sampler.
sumAll(&(pval_)[0],&pval_sum[0],pval_.size());
198 for (
unsigned p = 0; p < order_.size(); p++ ) {
199 dev += coeff_[p] * std::pow(pval_sum[p],1.0/order_[p]);
206 std::vector<Real> pval_sum(pval_.size());
207 sampler.
sumAll(&(pval_)[0],&pval_sum[0],pval_.size());
208 Teuchos::RCP<Vector<Real> > pg;
209 for (
unsigned p = 0; p < order_.size(); p++ ) {
210 if ( pval_sum[p] > 0.0 ) {
211 pg = (pg_[p])->clone();
212 sampler.
sumAll(*(pg_[p]),*pg);
213 g.
axpy(coeff_[p]/std::pow(pval_sum[p],1.0-1.0/order_[p]),*pg);
219 std::vector<Real> pval_sum(pval_.size());
220 sampler.
sumAll(&(pval_)[0],&pval_sum[0],pval_.size());
221 std::vector<Real> pgv_sum(pgv_.size());
222 sampler.
sumAll(&(pgv_)[0],&pgv_sum[0],pgv_.size());
224 Teuchos::RCP<Vector<Real> > pg, pg0, phv;
225 for (
unsigned p = 0; p < order_.size(); p++ ) {
226 if ( pval_sum[p] > 0.0 ) {
227 pg = (pg_[p])->clone();
228 sampler.
sumAll(*(pg_[p]),*pg);
229 pg0 = (pg0_[p])->clone();
230 sampler.
sumAll(*(pg0_[p]),*pg0);
231 phv = (phv_[p])->clone();
232 sampler.
sumAll(*(phv_[p]),*phv);
233 c = coeff_[p]*(pgv_sum[p]/std::pow(pval_sum[p],2.0-1.0/order_[p]));
235 c = coeff_[p]/std::pow(pval_sum[p],1.0-1.0/order_[p]);
void reset(Teuchos::RCP< Vector< Real > > &x0, 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...
MeanDeviationFromTarget(std::vector< Real > &target, std::vector< Real > &order, std::vector< Real > &coeff, Teuchos::RCP< PositiveFunction< Real > > &pf)
std::vector< Real > target_
void getGradient(Vector< Real > &g, SampleGenerator< Real > &sampler)
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
void update(const Real val, const Vector< Real > &g, const Real gv, const Vector< Real > &hv, const Real weight)
Defines the linear algebra or vector space interface.
void sumAll(Real *input, Real *output, int dim) const
std::vector< Real > coeff_
MeanDeviationFromTarget(Real target, Real order, Real coeff, Teuchos::RCP< PositiveFunction< Real > > &pf)
std::vector< Real > pval_
Real getValue(SampleGenerator< Real > &sampler)
std::vector< Teuchos::RCP< Vector< Real > > > phv_
void update(const Real val, const Real weight)
Teuchos::RCP< PositiveFunction< Real > > positiveFunction_
void reset(Teuchos::RCP< Vector< Real > > &x0, const Vector< Real > &x, Teuchos::RCP< Vector< Real > > &v0, const Vector< Real > &v)
virtual void reset(Teuchos::RCP< Vector< Real > > &x0, const Vector< Real > &x)
std::vector< Real > order_
std::vector< Teuchos::RCP< Vector< Real > > > pg_
std::vector< Teuchos::RCP< Vector< Real > > > pg0_
void getHessVec(Vector< Real > &hv, SampleGenerator< Real > &sampler)
void update(const Real val, const Vector< Real > &g, const Real weight)