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