ROL
ROL_Bundle_AS.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 #ifndef ROL_BUNDLE_AS_H
45 #define ROL_BUNDLE_AS_H
46 
47 #include "ROL_Bundle.hpp"
48 
53 namespace ROL {
54 
55 template<class Real>
56 class Bundle_AS : public Bundle<Real> {
57 /***********************************************************************************************/
58 /***************** BUNDLE STORAGE **************************************************************/
59 /***********************************************************************************************/
60 private:
61 
62  ROL::Ptr<Vector<Real> > tG_;
63  ROL::Ptr<Vector<Real> > eG_;
64  ROL::Ptr<Vector<Real> > yG_;
65  ROL::Ptr<Vector<Real> > gx_;
66  ROL::Ptr<Vector<Real> > ge_;
67 
68  std::set<unsigned> workingSet_;
69  std::set<unsigned> nworkingSet_;
70 
72 
73 /***********************************************************************************************/
74 /***************** BUNDLE MODIFICATION AND ACCESS ROUTINES *************************************/
75 /***********************************************************************************************/
76 public:
77  Bundle_AS(const unsigned maxSize = 10,
78  const Real coeff = 0.0,
79  const Real omega = 2.0,
80  const unsigned remSize = 2)
81  : Bundle<Real>(maxSize,coeff,omega,remSize), isInitialized_(false) {}
82 
83  void initialize(const Vector<Real> &g) {
85  if ( !isInitialized_ ) {
86  tG_ = g.clone();
87  yG_ = g.clone();
88  eG_ = g.clone();
89  gx_ = g.clone();
90  ge_ = g.clone();
91  isInitialized_ = true;
92  }
93  }
94 
95  unsigned solveDual(const Real t, const unsigned maxit = 1000, const Real tol = 1.e-8) {
96  unsigned iter = 0;
97  if (Bundle<Real>::size() == 1) {
98  iter = Bundle<Real>::solveDual_dim1(t,maxit,tol);
99  }
100  else if (Bundle<Real>::size() == 2) {
101  iter = Bundle<Real>::solveDual_dim2(t,maxit,tol);
102  }
103  else {
104  iter = solveDual_arbitrary(t,maxit,tol);
105  }
106  return iter;
107  }
108 
109 /***********************************************************************************************/
110 /***************** DUAL CUTTING PLANE PROBLEM ROUTINES *****************************************/
111 /***********************************************************************************************/
112 private:
113  void initializeDualSolver(void) {
114  Real sum(0), err(0), tmp(0), y(0), zero(0);
115  for (unsigned i = 0; i < Bundle<Real>::size(); ++i) {
116  // Compute sum of dualVariables_ using Kahan's compensated sum
117  //sum += Bundle<Real>::getDualVariable(i);
118  y = Bundle<Real>::getDualVariable(i) - err;
119  tmp = sum + y;
120  err = (tmp - sum) - y;
121  sum = tmp;
122  }
123  for (unsigned i = 0; i < Bundle<Real>::size(); ++i) {
124  tmp = Bundle<Real>::getDualVariable(i)/sum;
126  }
127  nworkingSet_.clear();
128  workingSet_.clear();
129  for (unsigned i = 0; i < Bundle<Real>::size(); ++i) {
130  if ( Bundle<Real>::getDualVariable(i) == zero ) {
131  workingSet_.insert(i);
132  }
133  else {
134  nworkingSet_.insert(i);
135  }
136  }
137  }
138 
139  void computeLagMult(std::vector<Real> &lam, const Real mu, const std::vector<Real> &g) const {
140  Real zero(0);
141  unsigned n = workingSet_.size();
142  if ( n > 0 ) {
143  lam.resize(n,zero);
144  typename std::set<unsigned>::iterator it = workingSet_.begin();
145  for (unsigned i = 0; i < n; ++i) {
146  lam[i] = g[*it] - mu; it++;
147  }
148  }
149  else {
150  lam.clear();
151  }
152  }
153 
154  bool isNonnegative(unsigned &ind, const std::vector<Real> &x) const {
155  bool flag = true;
156  unsigned n = workingSet_.size(); ind = Bundle<Real>::size();
157  if ( n > 0 ) {
158  Real min = ROL_OVERFLOW<Real>();
159  typename std::set<unsigned>::iterator it = workingSet_.begin();
160  for (unsigned i = 0; i < n; ++i) {
161  if ( x[i] < min ) {
162  ind = *it;
163  min = x[i];
164  }
165  it++;
166  }
167  flag = ((min < -ROL_EPSILON<Real>()) ? false : true);
168  }
169  return flag;
170  }
171 
172  Real computeStepSize(unsigned &ind, const std::vector<Real> &x, const std::vector<Real> &p) const {
173  Real alpha(1), tmp(0), zero(0); ind = Bundle<Real>::size();
174  typename std::set<unsigned>::iterator it;
175  for (it = nworkingSet_.begin(); it != nworkingSet_.end(); it++) {
176  if ( p[*it] < -ROL_EPSILON<Real>() ) {
177  tmp = -x[*it]/p[*it];
178  if ( alpha >= tmp ) {
179  alpha = tmp;
180  ind = *it;
181  }
182  }
183  }
184  return std::max(zero,alpha);
185  }
186 
187  unsigned solveEQPsubproblem(std::vector<Real> &s, Real &mu,
188  const std::vector<Real> &g, const Real tol) const {
189  // Build reduced QP information
190  Real zero(0);
191  unsigned n = nworkingSet_.size(), cnt = 0;
192  mu = zero;
193  s.assign(Bundle<Real>::size(),zero);
194  if ( n > 0 ) {
195  std::vector<Real> gk(n,zero);
196  typename std::set<unsigned>::iterator it = nworkingSet_.begin();
197  for (unsigned i = 0; i < n; ++i) {
198  gk[i] = g[*it]; it++;
199  }
200  std::vector<Real> sk(n,zero);
201  cnt = projectedCG(sk,mu,gk,tol);
202  it = nworkingSet_.begin();
203  for (unsigned i = 0; i < n; ++i) {
204  s[*it] = sk[i]; it++;
205  }
206  }
207  return cnt;
208  }
209 
210  void applyPreconditioner(std::vector<Real> &Px, const std::vector<Real> &x) const {
211  Real zero(0);
212  int type = 0;
213  std::vector<Real> tmp(Px.size(),zero);
214  switch (type) {
215  case 0: applyPreconditioner_Identity(tmp,x); break;
216  case 1: applyPreconditioner_Jacobi(tmp,x); break;
217  case 2: applyPreconditioner_SymGS(tmp,x); break;
218  }
220  }
221 
222  void applyG(std::vector<Real> &Gx, const std::vector<Real> &x) const {
223  int type = 0;
224  switch (type) {
225  case 0: applyG_Identity(Gx,x); break;
226  case 1: applyG_Jacobi(Gx,x); break;
227  case 2: applyG_SymGS(Gx,x); break;
228  }
229  }
230 
231  void applyPreconditioner_Identity(std::vector<Real> &Px, const std::vector<Real> &x) const {
232  unsigned dim = nworkingSet_.size();
233  Real sum(0), err(0), tmp(0), y(0);
234  for (unsigned i = 0; i < dim; ++i) {
235  // Compute sum of x using Kahan's compensated sum
236  //sum += x[i];
237  y = x[i] - err;
238  tmp = sum + y;
239  err = (tmp - sum) - y;
240  sum = tmp;
241  }
242  sum /= (Real)dim;
243  for (unsigned i = 0; i < dim; ++i) {
244  Px[i] = x[i] - sum;
245  }
246  }
247 
248  void applyG_Identity(std::vector<Real> &Gx, const std::vector<Real> &x) const {
249  Gx.assign(x.begin(),x.end());
250  }
251 
252  void applyPreconditioner_Jacobi(std::vector<Real> &Px, const std::vector<Real> &x) const {
253  unsigned dim = nworkingSet_.size();
254  Real eHe(0), sum(0), one(1), zero(0);
255  Real errX(0), tmpX(0), yX(0), errE(0), tmpE(0), yE(0);
256  std::vector<Real> gg(dim,zero);
257  typename std::set<unsigned>::iterator it = nworkingSet_.begin();
258  for (unsigned i = 0; i < dim; ++i) {
259  gg[i] = one/std::abs(Bundle<Real>::GiGj(*it,*it)); it++;
260  // Compute sum of inv(D)x using Kahan's aggregated sum
261  //sum += x[i]*gg[i];
262  yX = x[i]*gg[i] - errX;
263  tmpX = sum + yX;
264  errX = (tmpX - sum) - yX;
265  sum = tmpX;
266  // Compute sum of inv(D)e using Kahan's aggregated sum
267  //eHe += gg[i];
268  yE = gg[i] - errE;
269  tmpE = eHe + yE;
270  errE = (tmpE - eHe) - yE;
271  eHe = tmpE;
272  }
273  sum /= eHe;
274  for (unsigned i = 0; i < dim; ++i) {
275  Px[i] = (x[i]-sum)*gg[i];
276  }
277  }
278 
279  void applyG_Jacobi(std::vector<Real> &Gx, const std::vector<Real> &x) const {
280  unsigned dim = nworkingSet_.size();
281  typename std::set<unsigned>::iterator it = nworkingSet_.begin();
282  for (unsigned i = 0; i < dim; ++i) {
283  Gx[i] = std::abs(Bundle<Real>::GiGj(*it,*it))*x[i]; it++;
284  }
285  }
286 
287  void applyPreconditioner_SymGS(std::vector<Real> &Px, const std::vector<Real> &x) const {
288  int dim = nworkingSet_.size();
289  //unsigned cnt = 0;
290  gx_->zero(); ge_->zero();
291  Real eHx(0), eHe(0), one(1);
292  // Forward substitution
293  std::vector<Real> x1(dim,0), e1(dim,0),gg(dim,0);
294  typename std::set<unsigned>::iterator it, jt;
295  it = nworkingSet_.begin();
296  for (int i = 0; i < dim; ++i) {
297  gx_->zero(); ge_->zero(); jt = nworkingSet_.begin();
298  for (int j = 0; j < i; ++j) {
299  Bundle<Real>::addGi(*jt,x1[j],*gx_);
300  Bundle<Real>::addGi(*jt,e1[j],*ge_);
301  jt++;
302  }
303  gg[i] = Bundle<Real>::GiGj(*it,*it);
304  x1[i] = (x[i] - Bundle<Real>::dotGi(*it,*gx_))/gg[i];
305  e1[i] = (one - Bundle<Real>::dotGi(*it,*ge_))/gg[i];
306  it++;
307  }
308  // Apply diagonal
309  for (int i = 0; i < dim; ++i) {
310  x1[i] *= gg[i];
311  e1[i] *= gg[i];
312  }
313  // Back substitution
314  std::vector<Real> Hx(dim,0), He(dim,0); it = nworkingSet_.end();
315  for (int i = dim-1; i >= 0; --i) {
316  it--;
317  gx_->zero(); ge_->zero(); jt = nworkingSet_.end();
318  for (int j = dim-1; j >= i+1; --j) {
319  jt--;
320  Bundle<Real>::addGi(*jt,Hx[j],*gx_);
321  Bundle<Real>::addGi(*jt,He[j],*ge_);
322  }
323  Hx[i] = (x1[i] - Bundle<Real>::dotGi(*it,*gx_))/gg[i];
324  He[i] = (e1[i] - Bundle<Real>::dotGi(*it,*ge_))/gg[i];
325  // Compute sums
326  eHx += Hx[i];
327  eHe += He[i];
328  }
329  // Accumulate preconditioned vector
330  for (int i = 0; i < dim; ++i) {
331  Px[i] = Hx[i] - (eHx/eHe)*He[i];
332  }
333  }
334 
335  void applyG_SymGS(std::vector<Real> &Gx, const std::vector<Real> &x) const {
336  unsigned dim = nworkingSet_.size();
337  typename std::set<unsigned>::iterator it = nworkingSet_.begin();
338  for (unsigned i = 0; i < dim; ++i) {
339  Gx[i] = std::abs(Bundle<Real>::GiGj(*it,*it))*x[i]; it++;
340  }
341  }
342 
343  void computeResidualUpdate(std::vector<Real> &r, std::vector<Real> &g) const {
344  unsigned n = g.size();
345  std::vector<Real> Gg(n,0);
346  Real y(0), ytmp(0), yprt(0), yerr(0);
347  applyPreconditioner(g,r);
348  applyG(Gg,g);
349  // Compute multiplier using Kahan's compensated sum
350  for (unsigned i = 0; i < n; ++i) {
351  //y += (r[i] - Gg[i]);
352  yprt = (r[i] - Gg[i]) - yerr;
353  ytmp = y + yprt;
354  yerr = (ytmp - y) - yprt;
355  y = ytmp;
356  }
357  y /= (Real)n;
358  for (unsigned i = 0; i < n; ++i) {
359  r[i] -= y;
360  }
361  applyPreconditioner(g,r);
362  }
363 
364  void applyFullMatrix(std::vector<Real> &Hx, const std::vector<Real> &x) const {
365  Real one(1);
366  gx_->zero(); eG_->zero();
367  for (unsigned i = 0; i < Bundle<Real>::size(); ++i) {
368  // Compute Gx using Kahan's compensated sum
369  //gx_->axpy(x[i],Bundle<Real>::subgradient(i));
370  yG_->set(Bundle<Real>::subgradient(i)); yG_->scale(x[i]); yG_->axpy(-one,*eG_);
371  tG_->set(*gx_); tG_->plus(*yG_);
372  eG_->set(*tG_); eG_->axpy(-one,*gx_); eG_->axpy(-one,*yG_);
373  gx_->set(*tG_);
374  }
375  for (unsigned i = 0; i < Bundle<Real>::size(); ++i) {
376  // Compute < g_i, Gx >
377  Hx[i] = Bundle<Real>::dotGi(i,*gx_);
378  }
379  }
380 
381  void applyMatrix(std::vector<Real> &Hx, const std::vector<Real> &x) const {
382  Real one(1);
383  gx_->zero(); eG_->zero();
384  unsigned n = nworkingSet_.size();
385  typename std::set<unsigned>::iterator it = nworkingSet_.begin();
386  for (unsigned i = 0; i < n; ++i) {
387  // Compute Gx using Kahan's compensated sum
388  //gx_->axpy(x[i],Bundle<Real>::subgradient(*it));
389  yG_->set(Bundle<Real>::subgradient(*it)); yG_->scale(x[i]); yG_->axpy(-one,*eG_);
390  tG_->set(*gx_); tG_->plus(*yG_);
391  eG_->set(*tG_); eG_->axpy(-one,*gx_); eG_->axpy(-one,*yG_);
392  gx_->set(*tG_);
393  it++;
394  }
395  it = nworkingSet_.begin();
396  for (unsigned i = 0; i < n; ++i) {
397  // Compute < g_i, Gx >
398  Hx[i] = Bundle<Real>::dotGi(*it,*gx_); it++;
399  }
400  }
401 
402  unsigned projectedCG(std::vector<Real> &x, Real &mu, const std::vector<Real> &b, const Real tol) const {
403  Real one(1), zero(0);
404  unsigned n = nworkingSet_.size();
405  std::vector<Real> r(n,0), r0(n,0), g(n,0), d(n,0), Ad(n,0);
406  // Compute residual Hx + g = g with x = 0
407  x.assign(n,0);
408  scale(r,one,b);
409  r0.assign(r.begin(),r.end());
410  // Precondition residual
412  Real rg = dot(r,g), rg0(0);
413  // Get search direction
414  scale(d,-one,g);
415  Real alpha(0), kappa(0), beta(0), TOL(1.e-2);
416  Real CGtol = std::min(tol,TOL*rg);
417  unsigned cnt = 0;
418  while (rg > CGtol && cnt < 2*n+1) {
419  applyMatrix(Ad,d);
420  kappa = dot(d,Ad);
421  alpha = rg/kappa;
422  axpy(alpha,d,x);
423  axpy(alpha,Ad,r);
424  axpy(alpha,Ad,r0);
426  rg0 = rg;
427  rg = dot(r,g);
428  beta = rg/rg0;
429  scale(d,beta);
430  axpy(-one,g,d);
431  cnt++;
432  }
433  // Compute multiplier for equality constraint using Kahan's compensated sum
434  mu = zero;
435  Real err(0), tmp(0), y(0);
436  for (unsigned i = 0; i < n; ++i) {
437  //mu += r0[i];
438  y = r0[i] - err;
439  tmp = mu + y;
440  err = (tmp - mu) - y;
441  mu = tmp;
442  }
443  mu /= static_cast<Real>(n);
444  // Return iteration count
445  return cnt;
446  }
447 
448  Real dot(const std::vector<Real> &x, const std::vector<Real> &y) const {
449  // Compute dot product of two vectors using Kahan's compensated sum
450  Real val(0), err(0), tmp(0), y0(0);
451  unsigned n = std::min(x.size(),y.size());
452  for (unsigned i = 0; i < n; ++i) {
453  //val += x[i]*y[i];
454  y0 = x[i]*y[i] - err;
455  tmp = val + y0;
456  err = (tmp - val) - y0;
457  val = tmp;
458  }
459  return val;
460  }
461 
462  Real norm(const std::vector<Real> &x) const {
463  return std::sqrt(dot(x,x));
464  }
465 
466  void axpy(const Real a, const std::vector<Real> &x, std::vector<Real> &y) const {
467  unsigned n = std::min(y.size(),x.size());
468  for (unsigned i = 0; i < n; ++i) {
469  y[i] += a*x[i];
470  }
471  }
472 
473  void scale(std::vector<Real> &x, const Real a) const {
474  for (unsigned i = 0; i < x.size(); ++i) {
475  x[i] *= a;
476  }
477  }
478 
479  void scale(std::vector<Real> &x, const Real a, const std::vector<Real> &y) const {
480  unsigned n = std::min(x.size(),y.size());
481  for (unsigned i = 0; i < n; ++i) {
482  x[i] = a*y[i];
483  }
484  }
485 
486  unsigned solveDual_arbitrary(const Real t, const unsigned maxit = 1000, const Real tol = 1.e-8) {
488  bool nonneg = false;
489  unsigned ind = 0, i = 0, CGiter = 0;
490  Real snorm(0), alpha(0), mu(0), one(1), zero(0);
491  std::vector<Real> s(Bundle<Real>::size(),0), Hs(Bundle<Real>::size(),0);
492  std::vector<Real> g(Bundle<Real>::size(),0), lam(Bundle<Real>::size()+1,0);
493  std::vector<Real> dualVariables(Bundle<Real>::size(),0);
494  for (unsigned j = 0; j < Bundle<Real>::size(); ++j) {
495  dualVariables[j] = Bundle<Real>::getDualVariable(j);
496  }
497  //Real val = Bundle<Real>::evaluateObjective(g,dualVariables,t);
498  Bundle<Real>::evaluateObjective(g,dualVariables,t);
499  for (i = 0; i < maxit; ++i) {
500  CGiter += solveEQPsubproblem(s,mu,g,tol);
501  snorm = norm(s);
502  if ( snorm < ROL_EPSILON<Real>() ) {
503  computeLagMult(lam,mu,g);
504  nonneg = isNonnegative(ind,lam);
505  if ( nonneg ) {
506  break;
507  }
508  else {
509  alpha = one;
510  if ( ind < Bundle<Real>::size() ) {
511  workingSet_.erase(ind);
512  nworkingSet_.insert(ind);
513  }
514  }
515  }
516  else {
517  alpha = computeStepSize(ind,dualVariables,s);
518  if ( alpha > zero ) {
519  axpy(alpha,s,dualVariables);
520  applyFullMatrix(Hs,s);
521  axpy(alpha,Hs,g);
522  }
523  if (ind < Bundle<Real>::size()) {
524  workingSet_.insert(ind);
525  nworkingSet_.erase(ind);
526  }
527  }
528  //std::cout << "iter = " << i << " snorm = " << snorm << " alpha = " << alpha << "\n";
529  }
530  //Real crit = computeCriticality(g,dualVariables);
531  //std::cout << "Criticality Measure: " << crit << "\n";
532  //std::cout << "dim = " << Bundle<Real>::size() << " iter = " << i << " CGiter = " << CGiter << " CONVERGED!\n";
533  for (unsigned j = 0; j < Bundle<Real>::size(); ++j) {
534  Bundle<Real>::setDualVariable(j,dualVariables[j]);
535  }
536  return i;
537  }
538 
539  /************************************************************************/
540  /********************** PROJECTION ONTO FEASIBLE SET ********************/
541  /************************************************************************/
542  void project(std::vector<Real> &x, const std::vector<Real> &v) const {
543  std::vector<Real> vsort(Bundle<Real>::size(),0);
544  vsort.assign(v.begin(),v.end());
545  std::sort(vsort.begin(),vsort.end());
546  Real sum(-1), lam(0), zero(0), one(1);
547  for (unsigned i = Bundle<Real>::size()-1; i > 0; i--) {
548  sum += vsort[i];
549  if ( sum >= static_cast<Real>(Bundle<Real>::size()-i)*vsort[i-1] ) {
550  lam = sum/static_cast<Real>(Bundle<Real>::size()-i);
551  break;
552  }
553  }
554  if (lam == zero) {
555  lam = (sum+vsort[0])/static_cast<Real>(Bundle<Real>::size());
556  }
557  for (unsigned i = 0; i < Bundle<Real>::size(); ++i) {
558  x[i] = std::max(zero,v[i] - lam);
559  }
560  }
561 
562  Real computeCriticality(const std::vector<Real> &g, const std::vector<Real> &sol) {
563  Real zero(0), one(1);
564  std::vector<Real> x(Bundle<Real>::size(),0), Px(Bundle<Real>::size(),0);
565  axpy(one,sol,x);
566  axpy(-one,g,x);
567  project(Px,x);
568  scale(x,zero);
569  axpy(one,sol,x);
570  axpy(-one,Px,x);
571  return norm(x);
572  }
573 }; // class Bundle_AS
574 
575 } // namespace ROL
576 
577 #endif
std::set< unsigned > nworkingSet_
const Real dotGi(const unsigned i, const Vector< Real > &x) const
Definition: ROL_Bundle.hpp:292
ROL::Ptr< Vector< Real > > gx_
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
void applyPreconditioner(std::vector< Real > &Px, const std::vector< Real > &x) const
void applyG_Jacobi(std::vector< Real > &Gx, const std::vector< Real > &x) const
Real dot(const std::vector< Real > &x, const std::vector< Real > &y) const
void applyPreconditioner_SymGS(std::vector< Real > &Px, const std::vector< Real > &x) const
void applyPreconditioner_Jacobi(std::vector< Real > &Px, const std::vector< Real > &x) const
ROL::Ptr< Vector< Real > > eG_
void axpy(const Real a, const std::vector< Real > &x, std::vector< Real > &y) const
void applyG_Identity(std::vector< Real > &Gx, const std::vector< Real > &x) const
void setDualVariable(const unsigned i, const Real val)
Definition: ROL_Bundle.hpp:184
std::set< unsigned > workingSet_
Real norm(const std::vector< Real > &x) const
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
unsigned projectedCG(std::vector< Real > &x, Real &mu, const std::vector< Real > &b, const Real tol) const
void scale(std::vector< Real > &x, const Real a) const
Real computeCriticality(const std::vector< Real > &g, const std::vector< Real > &sol)
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Real computeStepSize(unsigned &ind, const std::vector< Real > &x, const std::vector< Real > &p) const
unsigned solveEQPsubproblem(std::vector< Real > &s, Real &mu, const std::vector< Real > &g, const Real tol) const
void computeLagMult(std::vector< Real > &lam, const Real mu, const std::vector< Real > &g) const
const Real alpha(const unsigned i) const
Definition: ROL_Bundle.hpp:201
void initialize(const Vector< Real > &g)
void computeResidualUpdate(std::vector< Real > &r, std::vector< Real > &g) const
ROL::Ptr< Vector< Real > > tG_
void scale(std::vector< Real > &x, const Real a, const std::vector< Real > &y) const
void applyFullMatrix(std::vector< Real > &Hx, const std::vector< Real > &x) const
bool isNonnegative(unsigned &ind, const std::vector< Real > &x) const
Real evaluateObjective(std::vector< Real > &g, const std::vector< Real > &x, const Real t) const
Definition: ROL_Bundle.hpp:300
const Real getDualVariable(const unsigned i) const
Definition: ROL_Bundle.hpp:180
virtual void initialize(const Vector< Real > &g)
Definition: ROL_Bundle.hpp:146
void initializeDualSolver(void)
Provides the interface for and implements an active set bundle.
void addGi(const unsigned i, const Real a, Vector< Real > &x) const
Definition: ROL_Bundle.hpp:296
unsigned solveDual_dim1(const Real t, const unsigned maxit=1000, const Real tol=1.e-8)
Definition: ROL_Bundle.hpp:327
ROL::Ptr< Vector< Real > > ge_
unsigned size(void) const
Definition: ROL_Bundle.hpp:205
void applyPreconditioner_Identity(std::vector< Real > &Px, const std::vector< Real > &x) const
void applyG(std::vector< Real > &Gx, const std::vector< Real > &x) const
const Real GiGj(const unsigned i, const unsigned j) const
Definition: ROL_Bundle.hpp:288
void applyMatrix(std::vector< Real > &Hx, const std::vector< Real > &x) const
unsigned solveDual_dim2(const Real t, const unsigned maxit=1000, const Real tol=1.e-8)
Definition: ROL_Bundle.hpp:333
void applyG_SymGS(std::vector< Real > &Gx, const std::vector< Real > &x) const
unsigned solveDual(const Real t, const unsigned maxit=1000, const Real tol=1.e-8)
ROL::Ptr< Vector< Real > > yG_
void project(std::vector< Real > &x, const std::vector< Real > &v) const
Provides the interface for and implements a bundle.
Definition: ROL_Bundle.hpp:62
unsigned solveDual_arbitrary(const Real t, const unsigned maxit=1000, const Real tol=1.e-8)
Bundle_AS(const unsigned maxSize=10, const Real coeff=0.0, const Real omega=2.0, const unsigned remSize=2)