54 #include "Teuchos_GlobalMPISession.hpp"
64 template <
class Real,
class Element=Real>
67 template <
class Real,
class Element=Real>
70 template <
class Real,
class Element=Real>
73 template <
class Real,
class Element=Real>
79 template <
class Real,
class Element>
87 ROL::Ptr<std::vector<Element> >
std_vec_;
88 mutable ROL::Ptr<OptDualStdVector<Real> >
dual_vec_;
101 (*std_vec_)[i] += (*xp)[i];
108 (*std_vec_)[i] *= alpha;
121 val += (*std_vec_)[i]*(*xp)[i];
128 val = std::sqrt(
dot(*
this) );
132 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
133 return ROL::makePtr<OptStdVector>( ROL::makePtr<std::vector<Element>>(
std_vec_->size()) );
136 ROL::Ptr<const std::vector<Element> >
getVector()
const {
144 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
146 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(
std_vec_->size(),0.0);
147 ROL::Ptr<V> e = ROL::makePtr<OptStdVector>( e_ptr );
155 dual_vec_ = ROL::makePtr<OptDualStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
163 template <
class Real,
class Element>
171 ROL::Ptr<std::vector<Element> >
std_vec_;
172 mutable ROL::Ptr<OptStdVector<Real> >
dual_vec_;
183 (*std_vec_)[i] += (*xp)[i];
190 (*std_vec_)[i] *= alpha;
200 val += (*std_vec_)[i]*(*xp)[i];
207 val = std::sqrt(
dot(*
this) );
211 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
212 return ROL::makePtr<OptDualStdVector>( ROL::makePtr<std::vector<Element>>(
std_vec_->size()) );
215 ROL::Ptr<const std::vector<Element> >
getVector()
const {
223 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
225 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(
std_vec_->size(), 0.0 );
226 ROL::Ptr<V> e = ROL::makePtr<OptDualStdVector>( e_ptr );
234 dual_vec_ = ROL::makePtr<OptStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
242 template <
class Real,
class Element>
262 (*std_vec_)[i] += (*xp)[i];
269 (*std_vec_)[i] *= alpha;
279 val += (*std_vec_)[i]*(*xp)[i];
286 val = std::sqrt(
dot(*
this) );
290 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
291 return ROL::makePtr<ConStdVector>( ROL::makePtr<std::vector<Element>>(
std_vec_->size()));
294 ROL::Ptr<const std::vector<Element> >
getVector()
const {
302 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
304 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(
std_vec_->size(), 0.0);
305 ROL::Ptr<V> e = ROL::makePtr<ConStdVector>( e_ptr );
313 dual_vec_ = ROL::makePtr<ConDualStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
321 template <
class Real,
class Element>
342 (*std_vec_)[i] += (*xp)[i];
349 (*std_vec_)[i] *= alpha;
359 val += (*std_vec_)[i]*(*xp)[i];
366 val = std::sqrt(
dot(*
this) );
370 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
371 return ROL::makePtr<ConDualStdVector>( ROL::makePtr<std::vector<Element>>(
std_vec_->size()));
374 ROL::Ptr<const std::vector<Element> >
getVector()
const {
382 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
384 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(
std_vec_->size(),0.0);
385 ROL::Ptr<V> e = ROL::makePtr<ConDualStdVector>(e_ptr);
393 dual_vec_ = ROL::makePtr<ConStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
402 int main(
int argc,
char *argv[]) {
404 typedef std::vector<RealT> vector;
409 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
412 int iprint = argc - 1;
413 ROL::Ptr<std::ostream> outStream;
416 outStream = ROL::makePtrFromRef(std::cout);
418 outStream = ROL::makePtrFromRef(bhs);
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);
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);
452 OptDualStdVector<RealT> g(g_ptr);
454 OptDualStdVector<RealT> gd(gd_ptr);
456 ConStdVector<RealT> vc(vc_ptr);
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;
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;
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";
477 ROL::Ptr<vector> v1_ptr = ROL::makePtr<vector>(
dim, 0.0);
478 ROL::Ptr<vector> v2_ptr = ROL::makePtr<vector>(nc, 0.0);
482 constr->solveAugmentedSystem(v1, v2, gd, vc, xtest, augtol);
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);
507 algo.run(x, g, vl, vc, *obj, *constr,
true, *outStream);
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";
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>()) ) {
531 catch (std::logic_error& err) {
532 *outStream << err.what() <<
"\n";
537 std::cout <<
"End Result: TEST FAILED\n";
539 std::cout <<
"End Result: TEST PASSED\n";
std::vector< Element > vector
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 .
std::vector< Element > vector
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
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.
Real norm() const
Returns where .
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.
ROL::Ptr< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
std::vector< Element > vector
Defines the linear algebra or vector space interface.
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.
std::vector< Element > 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.
Real norm() const
Returns where .
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
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 .
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()
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 .