ROL
dual-spaces/simple-eq-constr-1/example_01.cpp
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 
16 #include "ROL_Solver.hpp"
17 #include "ROL_Stream.hpp"
18 #include "Teuchos_GlobalMPISession.hpp"
19 
20 #include <iostream>
21 //#include <fenv.h>
22 
23 typedef double RealT;
24 
25 /*** Declare four vector spaces. ***/
26 
27 // Forward declarations:
28 
29 template <class Real, class Element=Real>
30 class OptStdVector; // Optimization space.
31 
32 template <class Real, class Element=Real>
33 class OptDualStdVector; // Dual optimization space.
34 
35 template <class Real, class Element=Real>
36 class ConStdVector; // Constraint space.
37 
38 template <class Real, class Element=Real>
39 class ConDualStdVector; // Dual constraint space.
40 
41 // Vector space definitions:
42 
43 // Optimization space.
44 template <class Real, class Element>
45 class OptStdVector : public ROL::Vector<Real> {
46 
47 typedef std::vector<Element> vector;
49 typedef typename vector::size_type uint;
50 
51 private:
52 ROL::Ptr<std::vector<Element> > std_vec_;
53 mutable ROL::Ptr<OptDualStdVector<Real> > dual_vec_;
54 
55 public:
56 
57 OptStdVector(const ROL::Ptr<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(ROL::nullPtr) {}
58 
59 void plus( const ROL::Vector<Real> &x ) {
60 
61 
62  ROL::Ptr<const vector> xp = dynamic_cast<const OptStdVector&>(x).getVector();
63 
64  uint dimension = std_vec_->size();
65  for (uint i=0; i<dimension; i++) {
66  (*std_vec_)[i] += (*xp)[i];
67  }
68 }
69 
70 void scale( const Real alpha ) {
71  uint dimension = std_vec_->size();
72  for (uint i=0; i<dimension; i++) {
73  (*std_vec_)[i] *= alpha;
74  }
75 }
76 
77 Real dot( const ROL::Vector<Real> &x ) const {
78 
79 
80 
81  ROL::Ptr<const vector> xp = dynamic_cast<const OptStdVector&>(x).getVector();
82  Real val = 0;
83 
84  uint dimension = std_vec_->size();
85  for (uint i=0; i<dimension; i++) {
86  val += (*std_vec_)[i]*(*xp)[i];
87  }
88  return val;
89 }
90 
91 Real norm() const {
92  Real val = 0;
93  val = std::sqrt( dot(*this) );
94  return val;
95 }
96 
97 ROL::Ptr<ROL::Vector<Real> > clone() const {
98  return ROL::makePtr<OptStdVector>( ROL::makePtr<std::vector<Element>>(std_vec_->size()) );
99 }
100 
101 ROL::Ptr<const std::vector<Element> > getVector() const {
102  return std_vec_;
103 }
104 
105 ROL::Ptr<std::vector<Element> > getVector() {
106  return std_vec_;
107 }
108 
109 ROL::Ptr<ROL::Vector<Real> > basis( const int i ) const {
110 
111  ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(std_vec_->size(),0.0);
112  ROL::Ptr<V> e = ROL::makePtr<OptStdVector>( e_ptr );
113  (*e_ptr)[i] = 1.0;
114  return e;
115 }
116 
117 int dimension() const {return static_cast<int>(std_vec_->size());}
118 
119 const ROL::Vector<Real> & dual() const {
120  dual_vec_ = ROL::makePtr<OptDualStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
121  return *dual_vec_;
122 }
123 
124 Real apply( const ROL::Vector<Real> &x ) const {
125  Real val = 0;
126 
127  ROL::Ptr<const vector> xvalptr = dynamic_cast<const OptDualStdVector<Real,Element>&>(x).getVector();
128  uint dimension = std_vec_->size();
129  for (uint i=0; i<dimension; i++) {
130  val += (*std_vec_)[i]*(*xvalptr)[i];
131  }
132  return val;
133 }
134 
135 }; // class OptStdVector
136 
137 
138 // Dual optimization space.
139 template <class Real, class Element>
140 class OptDualStdVector : public ROL::Vector<Real> {
141 
142 typedef std::vector<Element> vector;
144 typedef typename vector::size_type uint;
145 
146 private:
147 ROL::Ptr<std::vector<Element> > std_vec_;
148 mutable ROL::Ptr<OptStdVector<Real> > dual_vec_;
149 
150 public:
151 
152 OptDualStdVector(const ROL::Ptr<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(ROL::nullPtr) {}
153 
154 void plus( const ROL::Vector<Real> &x ) {
155 
156  ROL::Ptr<const vector> xp = dynamic_cast<const OptDualStdVector&>(x).getVector();
157  uint dimension = std_vec_->size();
158  for (uint i=0; i<dimension; i++) {
159  (*std_vec_)[i] += (*xp)[i];
160  }
161 }
162 
163 void scale( const Real alpha ) {
164  uint dimension = std_vec_->size();
165  for (uint i=0; i<dimension; i++) {
166  (*std_vec_)[i] *= alpha;
167  }
168 }
169 
170 Real dot( const ROL::Vector<Real> &x ) const {
171 
172  ROL::Ptr<const vector> xp = dynamic_cast<const OptDualStdVector&>(x).getVector();
173  Real val = 0;
174  uint dimension = std_vec_->size();
175  for (uint i=0; i<dimension; i++) {
176  val += (*std_vec_)[i]*(*xp)[i];
177  }
178  return val;
179 }
180 
181 Real norm() const {
182  Real val = 0;
183  val = std::sqrt( dot(*this) );
184  return val;
185 }
186 
187 ROL::Ptr<ROL::Vector<Real> > clone() const {
188  return ROL::makePtr<OptDualStdVector>( ROL::makePtr<std::vector<Element>>(std_vec_->size()) );
189 }
190 
191 ROL::Ptr<const std::vector<Element> > getVector() const {
192  return std_vec_;
193 }
194 
195 ROL::Ptr<std::vector<Element> > getVector() {
196  return std_vec_;
197 }
198 
199 ROL::Ptr<ROL::Vector<Real> > basis( const int i ) const {
200 
201  ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>( std_vec_->size(), 0.0 );
202  ROL::Ptr<V> e = ROL::makePtr<OptDualStdVector>( e_ptr );
203  (*e_ptr)[i] = 1.0;
204  return e;
205 }
206 
207 int dimension() const {return static_cast<int>(std_vec_->size());}
208 
209 const ROL::Vector<Real> & dual() const {
210  dual_vec_ = ROL::makePtr<OptStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
211  return *dual_vec_;
212 }
213 
214 Real apply( const ROL::Vector<Real> &x ) const {
215  Real val = 0;
216 
217  ROL::Ptr<const vector> xvalptr = dynamic_cast<const OptStdVector<Real,Element>&>(x).getVector();
218  uint dimension = std_vec_->size();
219  for (uint i=0; i<dimension; i++) {
220  val += (*std_vec_)[i]*(*xvalptr)[i];
221  }
222  return val;
223 }
224 
225 }; // class OptDualStdVector
226 
227 
228 // Constraint space.
229 template <class Real, class Element>
230 class ConStdVector : public ROL::Vector<Real> {
231 
232 typedef std::vector<Element> vector;
234 typedef typename vector::size_type uint;
235 
236 private:
237 ROL::Ptr<std::vector<Element> > std_vec_;
238 mutable ROL::Ptr<ConDualStdVector<Real> > dual_vec_;
239 
240 public:
241 
242 ConStdVector(const ROL::Ptr<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(ROL::nullPtr) {}
243 
244 void plus( const ROL::Vector<Real> &x ) {
245 
246  ROL::Ptr<const vector> xp = dynamic_cast<const ConStdVector&>(x).getVector();
247  uint dimension = std_vec_->size();
248  for (uint i=0; i<dimension; i++) {
249  (*std_vec_)[i] += (*xp)[i];
250  }
251 }
252 
253 void scale( const Real alpha ) {
254  uint dimension = std_vec_->size();
255  for (uint i=0; i<dimension; i++) {
256  (*std_vec_)[i] *= alpha;
257  }
258 }
259 
260 Real dot( const ROL::Vector<Real> &x ) const {
261 
262  ROL::Ptr<const vector> xp = dynamic_cast<const ConStdVector&>(x).getVector();
263  Real val = 0;
264  uint dimension = std_vec_->size();
265  for (uint i=0; i<dimension; i++) {
266  val += (*std_vec_)[i]*(*xp)[i];
267  }
268  return val;
269 }
270 
271 Real norm() const {
272  Real val = 0;
273  val = std::sqrt( dot(*this) );
274  return val;
275 }
276 
277 ROL::Ptr<ROL::Vector<Real> > clone() const {
278  return ROL::makePtr<ConStdVector>( ROL::makePtr<std::vector<Element>>(std_vec_->size()));
279 }
280 
281 ROL::Ptr<const std::vector<Element> > getVector() const {
282  return std_vec_;
283 }
284 
285 ROL::Ptr<std::vector<Element> > getVector() {
286  return std_vec_;
287 }
288 
289 ROL::Ptr<ROL::Vector<Real> > basis( const int i ) const {
290 
291  ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(std_vec_->size(), 0.0);
292  ROL::Ptr<V> e = ROL::makePtr<ConStdVector>( e_ptr );
293  (*e_ptr)[i] = 1.0;
294  return e;
295 }
296 
297 int dimension() const {return static_cast<int>(std_vec_->size());}
298 
299 const ROL::Vector<Real> & dual() const {
300  dual_vec_ = ROL::makePtr<ConDualStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
301  return *dual_vec_;
302 }
303 
304 Real apply( const ROL::Vector<Real> &x ) const {
305  Real val = 0;
306 
307  ROL::Ptr<const vector> xvalptr = dynamic_cast<const ConDualStdVector<Real,Element>&>(x).getVector();
308  uint dimension = std_vec_->size();
309  for (uint i=0; i<dimension; i++) {
310  val += (*std_vec_)[i]*(*xvalptr)[i];
311  }
312  return val;
313 }
314 
315 }; // class ConStdVector
316 
317 
318 // Dual constraint space.
319 template <class Real, class Element>
320 class ConDualStdVector : public ROL::Vector<Real> {
321 
322  typedef std::vector<Element> vector;
324  typedef typename vector::size_type uint;
325 
326 private:
327 
328 ROL::Ptr<std::vector<Element> > std_vec_;
329 mutable ROL::Ptr<ConStdVector<Real> > dual_vec_;
330 
331 public:
332 
333 ConDualStdVector(const ROL::Ptr<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(ROL::nullPtr) {}
334 
335 void plus( const ROL::Vector<Real> &x ) {
336 
337  ROL::Ptr<const vector> xp = dynamic_cast<const ConDualStdVector&>(x).getVector();
338  uint dimension = std_vec_->size();
339  for (uint i=0; i<dimension; i++) {
340  (*std_vec_)[i] += (*xp)[i];
341  }
342 }
343 
344 void scale( const Real alpha ) {
345  uint dimension = std_vec_->size();
346  for (uint i=0; i<dimension; i++) {
347  (*std_vec_)[i] *= alpha;
348  }
349 }
350 
351 Real dot( const ROL::Vector<Real> &x ) const {
352 
353  ROL::Ptr<const vector> xp = dynamic_cast<const ConDualStdVector&>(x).getVector();
354  Real val = 0;
355  uint dimension = std_vec_->size();
356  for (uint i=0; i<dimension; i++) {
357  val += (*std_vec_)[i]*(*xp)[i];
358  }
359  return val;
360 }
361 
362 Real norm() const {
363  Real val = 0;
364  val = std::sqrt( dot(*this) );
365  return val;
366 }
367 
368 ROL::Ptr<ROL::Vector<Real> > clone() const {
369  return ROL::makePtr<ConDualStdVector>( ROL::makePtr<std::vector<Element>>(std_vec_->size()));
370 }
371 
372 ROL::Ptr<const std::vector<Element> > getVector() const {
373  return std_vec_;
374 }
375 
376 ROL::Ptr<std::vector<Element> > getVector() {
377  return std_vec_;
378 }
379 
380 ROL::Ptr<ROL::Vector<Real> > basis( const int i ) const {
381 
382  ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(std_vec_->size(),0.0);
383  ROL::Ptr<V> e = ROL::makePtr<ConDualStdVector>(e_ptr);
384  (*e_ptr)[i] = 1.0;
385  return e;
386 }
387 
388 int dimension() const {return static_cast<int>(std_vec_->size());}
389 
390 const ROL::Vector<Real> & dual() const {
391  dual_vec_ = ROL::makePtr<ConStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
392  return *dual_vec_;
393 }
394 
395 Real apply( const ROL::Vector<Real> &x ) const {
396  Real val = 0;
397 
398  ROL::Ptr<const vector> xvalptr = dynamic_cast<const ConStdVector<Real,Element>&>(x).getVector();
399  uint dimension = std_vec_->size();
400  for (uint i=0; i<dimension; i++) {
401  val += (*std_vec_)[i]*(*xvalptr)[i];
402  }
403  return val;
404 }
405 
406 }; // class ConDualStdVector
407 
408 /*** End of declaration of four vector spaces. ***/
409 
410 
411 int main(int argc, char *argv[]) {
412  //feenableexcept(FE_DIVBYZERO | FE_INVALID | FE_OVERFLOW);
413 
414  typedef std::vector<RealT> vector;
415  typedef vector::size_type uint;
416 
417  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
418 
419  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
420  int iprint = argc - 1;
421  ROL::Ptr<std::ostream> outStream;
422  ROL::nullstream bhs; // outputs nothing
423  if (iprint > 0)
424  outStream = ROL::makePtrFromRef(std::cout);
425  else
426  outStream = ROL::makePtrFromRef(bhs);
427 
428  int errorFlag = 0;
429 
430  // *** Example body.
431 
432  try {
433 
434  uint dim = 5;
435  uint nc = 3;
436  ROL::Ptr<ROL::Objective<RealT> > obj;
437  ROL::Ptr<ROL::Constraint<RealT> > constr;
438  ROL::Ptr<vector> x_ptr = ROL::makePtr<vector>(dim, 0.0);
439  ROL::Ptr<vector> sol_ptr = ROL::makePtr<vector>(dim, 0.0);
440  ROL::Ptr<OptStdVector<RealT>> x = ROL::makePtr<OptStdVector<RealT>>(x_ptr); // Iteration vector.
441  ROL::Ptr<OptStdVector<RealT>> sol = ROL::makePtr<OptStdVector<RealT>>(sol_ptr); // Reference solution vector.
442 
443  // Retrieve objective, constraint, iteration vector, solution vector.
445  obj = SEC.getObjective();
446  constr = SEC.getEqualityConstraint();
447  x->set(*SEC.getInitialGuess());
448  sol->set(*SEC.getSolution());
449 
450  // Run derivative checks, etc.
451  RealT left = -1e0, right = 1e0;
452  ROL::Ptr<vector> xtest_ptr = ROL::makePtr<vector>(dim, 0.0);
453  ROL::Ptr<vector> g_ptr = ROL::makePtr<vector>(dim, 0.0);
454  ROL::Ptr<vector> d_ptr = ROL::makePtr<vector>(dim, 0.0);
455  ROL::Ptr<vector> gd_ptr = ROL::makePtr<vector>(dim, 0.0);
456  ROL::Ptr<vector> v_ptr = ROL::makePtr<vector>(dim, 0.0);
457  ROL::Ptr<vector> vc_ptr = ROL::makePtr<vector>(nc, 0.0);
458  ROL::Ptr<vector> vl_ptr = ROL::makePtr<vector>(nc, 0.0);
459  ROL::Ptr<OptStdVector<RealT>> xtest = ROL::makePtr<OptStdVector<RealT>>(xtest_ptr);
460  ROL::Ptr<OptDualStdVector<RealT>> g = ROL::makePtr<OptDualStdVector<RealT>>(g_ptr);
461  ROL::Ptr<OptStdVector<RealT>> d = ROL::makePtr<OptStdVector<RealT>>(d_ptr);
462  ROL::Ptr<OptDualStdVector<RealT>> gd = ROL::makePtr<OptDualStdVector<RealT>>(gd_ptr);
463  ROL::Ptr<OptStdVector<RealT>> v = ROL::makePtr<OptStdVector<RealT>>(v_ptr);
464  ROL::Ptr<ConStdVector<RealT>> vc = ROL::makePtr<ConStdVector<RealT>>(vc_ptr);
465  ROL::Ptr<ConDualStdVector<RealT>> vl = ROL::makePtr<ConDualStdVector<RealT>>(vl_ptr);
466  // set xtest, d, v
467  for (uint i=0; i<dim; i++) {
468  (*xtest_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
469  (*d_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
470  (*gd_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
471  (*v_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
472  }
473  // set vc, vl
474  for (uint i=0; i<nc; i++) {
475  (*vc_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
476  (*vl_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
477  }
478  obj->checkGradient(*xtest, *g, *d, true, *outStream); *outStream << std::endl;
479  obj->checkHessVec(*xtest, *g, *v, true, *outStream); *outStream << std::endl;
480  obj->checkHessSym(*xtest, *g, *d, *v, true, *outStream); *outStream << std::endl;
481  constr->checkApplyJacobian(*xtest, *v, *vc, true, *outStream); *outStream << std::endl;
482  constr->checkApplyAdjointJacobian(*xtest, *vl, *vc, *g, true, *outStream); *outStream << std::endl;
483  constr->checkApplyAdjointHessian(*xtest, *vl, *d, *g, true, *outStream); *outStream << std::endl;
484 
485  ROL::Ptr<vector> v1_ptr = ROL::makePtr<vector>(dim, 0.0);
486  ROL::Ptr<vector> v2_ptr = ROL::makePtr<vector>(nc, 0.0);
487  ROL::Ptr<OptStdVector<RealT>> v1 = ROL::makePtr<OptStdVector<RealT>>(v1_ptr);
488  ROL::Ptr<ConDualStdVector<RealT>> v2 = ROL::makePtr<ConDualStdVector<RealT>>(v2_ptr);
489  RealT augtol = 1e-8;
490  constr->solveAugmentedSystem(*v1, *v2, *gd, *vc, *xtest, augtol);
491 
492  // Define problem.
493  vl->zero(); vc->zero(); g->zero();
494  ROL::Ptr<ROL::Problem<RealT>>
495  problem = ROL::makePtr<ROL::Problem<RealT>>(obj,x,g);
496  problem->addConstraint("Equality",constr,ROL::staticPtrCast<ROL::Vector<RealT>>(vl),ROL::staticPtrCast<ROL::Vector<RealT>>(vc));
497  problem->finalize(false,true,*outStream);
498 
499  // Define algorithm.
500  ROL::ParameterList parlist;
501  std::string stepname = "Augmented Lagrangian";
502  parlist.sublist("General").set("Output Level",1);
503  parlist.sublist("Step").set("Type",stepname);
504  parlist.sublist("Step").sublist(stepname).sublist("Optimality System Solver").set("Nominal Relative Tolerance",1e-4);
505  parlist.sublist("Step").sublist(stepname).sublist("Optimality System Solver").set("Fix Tolerance",true);
506  parlist.sublist("Step").sublist(stepname).sublist("Tangential Subproblem Solver").set("Iteration Limit",20);
507  parlist.sublist("Step").sublist(stepname).sublist("Tangential Subproblem Solver").set("Relative Tolerance",1e-2);
508  parlist.sublist("Step").sublist(stepname).set("Output Level",0);
509  parlist.sublist("Status Test").set("Gradient Tolerance",1.e-12);
510  parlist.sublist("Status Test").set("Constraint Tolerance",1.e-12);
511  parlist.sublist("Status Test").set("Step Tolerance",1.e-18);
512  parlist.sublist("Status Test").set("Iteration Limit",100);
513  ROL::Ptr<ROL::Solver<RealT>> solver
514  = ROL::makePtr<ROL::Solver<RealT>>(problem,parlist);
515 
516  // Run Algorithm
517  //(*x_ptr)[0] = 3.0; (*x_ptr)[1] = 2.0; (*x_ptr)[2] = 2.0; (*x_ptr)[3] = 1.0; (*x_ptr)[4] = 1.0;
518  //(*x_ptr)[0] = -5.0; (*x_ptr)[1] = -5.0; (*x_ptr)[2] = -5.0; (*x_ptr)[3] = -6.0; (*x_ptr)[4] = -6.0;
519  solver->solve(*outStream);
520 
521  // Compute Error
522  *outStream << "\nReference solution x_r =\n";
523  *outStream << std::scientific << " " << (*sol_ptr)[0] << "\n";
524  *outStream << std::scientific << " " << (*sol_ptr)[1] << "\n";
525  *outStream << std::scientific << " " << (*sol_ptr)[2] << "\n";
526  *outStream << std::scientific << " " << (*sol_ptr)[3] << "\n";
527  *outStream << std::scientific << " " << (*sol_ptr)[4] << "\n";
528  *outStream << "\nOptimal solution x =\n";
529  *outStream << std::scientific << " " << (*x_ptr)[0] << "\n";
530  *outStream << std::scientific << " " << (*x_ptr)[1] << "\n";
531  *outStream << std::scientific << " " << (*x_ptr)[2] << "\n";
532  *outStream << std::scientific << " " << (*x_ptr)[3] << "\n";
533  *outStream << std::scientific << " " << (*x_ptr)[4] << "\n";
534  x->axpy(-1.0, *sol);
535  RealT abserr = x->norm();
536  RealT relerr = abserr/sol->norm();
537  *outStream << std::scientific << "\n Absolute Error: " << abserr;
538  *outStream << std::scientific << "\n Relative Error: " << relerr << "\n";
539  if ( relerr > sqrt(ROL::ROL_EPSILON<RealT>()) ) {
540  errorFlag += 1;
541  }
542  }
543  catch (std::logic_error& err) {
544  *outStream << err.what() << "\n";
545  errorFlag = -1000;
546  }; // end try
547 
548  if (errorFlag != 0)
549  std::cout << "End Result: TEST FAILED\n";
550  else
551  std::cout << "End Result: TEST PASSED\n";
552 
553  return 0;
554 
555 }
556 
ROL::Ptr< std::vector< Element > > std_vec_
typename PV< Real >::size_type size_type
void scale(const Real alpha)
Compute where .
Real dot(const ROL::Vector< Real > &x) const
Compute where .
ROL::Ptr< const std::vector< Element > > getVector() const
ROL::Ptr< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
ROL::Ptr< std::vector< Element > > getVector()
OptDualStdVector(const ROL::Ptr< std::vector< Element > > &std_vec)
ROL::Ptr< std::vector< Element > > std_vec_
ROL::Ptr< std::vector< Element > > getVector()
int dimension() const
Return dimension of the vector space.
ROL::Ptr< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:46
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
Real apply(const ROL::Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
ROL::Ptr< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
void plus(const ROL::Vector< Real > &x)
Compute , where .
ROL::Ptr< OptStdVector< Real > > dual_vec_
ROL::Ptr< std::vector< Element > > std_vec_
ROL::Ptr< std::vector< Element > > getVector()
basic_nullstream< char, std::char_traits< char >> nullstream
Definition: ROL_Stream.hpp:36
ROL::Ptr< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
ROL::Ptr< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
void scale(const Real alpha)
Compute where .
ConDualStdVector(const ROL::Ptr< std::vector< Element > > &std_vec)
void scale(const Real alpha)
Compute where .
ROL::Ptr< const std::vector< Element > > getVector() const
Real apply(const ROL::Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
Real dot(const ROL::Vector< Real > &x) const
Compute where .
int dimension() const
Return dimension of the vector space.
int dimension() const
Return dimension of the vector space.
ROL::Ptr< ConStdVector< Real > > dual_vec_
ROL::Ptr< OptDualStdVector< Real > > dual_vec_
Real apply(const ROL::Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
ConStdVector(const ROL::Ptr< std::vector< Element > > &std_vec)
ROL::Ptr< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Real apply(const ROL::Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
Contains definitions for the equality constrained NLP from Nocedal/Wright, 2nd edition, page 574, example 18.2; note the typo in reversing the initial guess and the solution.
int dimension() const
Return dimension of the vector space.
ROL::Ptr< std::vector< Element > > std_vec_
int main(int argc, char *argv[])
Real norm() const
Returns where .
ROL::Ptr< const std::vector< Element > > getVector() const
Real dot(const ROL::Vector< Real > &x) const
Compute where .
Real dot(const ROL::Vector< Real > &x) const
Compute where .
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:175
ROL::Ptr< const std::vector< Element > > getVector() const
ROL::Ptr< ConDualStdVector< Real > > dual_vec_
ROL::Ptr< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
void plus(const ROL::Vector< Real > &x)
Compute , where .
ROL::Ptr< std::vector< Element > > getVector()
constexpr auto dim
void plus(const ROL::Vector< Real > &x)
Compute , where .
void plus(const ROL::Vector< Real > &x)
Compute , where .
Real norm() const
Returns where .
OptStdVector(const ROL::Ptr< std::vector< Element > > &std_vec)
ROL::Ptr< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
void scale(const Real alpha)
Compute where .