20 #include "Teuchos_GlobalMPISession.hpp"
30 template <
class Real,
class Element=Real>
33 template <
class Real,
class Element=Real>
40 template <
class Real,
class Element>
49 ROL::Ptr<std::vector<Element> >
std_vec_;
50 mutable ROL::Ptr<OptDualStdVector<Real> >
dual_vec_;
61 (*std_vec_)[i] += (*xvalptr)[i];
65 void scale(
const Real alpha ) {
68 (*std_vec_)[i] *= alpha;
78 val += (*std_vec_)[i]*(*xvalptr)[i];
85 val = std::sqrt(
dot(*
this) );
89 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
90 return ROL::makePtr<OptStdVector>( ROL::makePtr<std::vector<Element>>(
std_vec_->size()) );
93 ROL::Ptr<const std::vector<Element> >
getVector()
const {
101 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
103 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(
std_vec_->size(),0.0);
106 ROL::Ptr<V> e = ROL::makePtr<OptStdVector>( e_ptr );
114 dual_vec_ = ROL::makePtr<OptDualStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
124 val += (*std_vec_)[i]*(*xvalptr)[i];
129 void applyUnary(
const ROL::Elementwise::UnaryFunction<Real> &f ) {
130 for(
auto& e : *
std_vec_ ) e = f.apply(e);
137 (*std_vec_)[i] = f.apply((*
std_vec_)[i],(*xvalptr)[i]);
141 Real
reduce(
const ROL::Elementwise::ReductionOp<Real> &r )
const {
142 Real result = r.initialValue();
154 template <
class Real,
class Element>
163 ROL::Ptr<std::vector<Element> >
std_vec_;
164 mutable ROL::Ptr<OptStdVector<Real> >
dual_vec_;
176 (*std_vec_)[i] += (*xvalptr)[i];
183 (*std_vec_)[i] *= alpha;
193 val += (*std_vec_)[i]*(*xvalptr)[i];
200 val = std::sqrt(
dot(*
this) );
204 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
205 return ROL::makePtr<OptDualStdVector>( ROL::makePtr<std::vector<Element>>(
std_vec_->size()) );
208 ROL::Ptr<const std::vector<Element> >
getVector()
const {
216 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
218 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(
std_vec_->size(),0.0);
221 ROL::Ptr<V> e = ROL::makePtr<OptDualStdVector>( e_ptr );
228 dual_vec_ = ROL::makePtr<OptStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
238 val += (*std_vec_)[i]*(*xvalptr)[i];
243 void applyUnary(
const ROL::Elementwise::UnaryFunction<Real> &f ) {
244 for(
auto& e : *
std_vec_ ) e = f.apply(e);
251 (*std_vec_)[i] = f.apply((*
std_vec_)[i],(*xvalptr)[i]);
255 Real
reduce(
const ROL::Elementwise::ReductionOp<Real> &r )
const {
256 Real result = r.initialValue();
274 int main(
int argc,
char *argv[]) {
276 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
279 int iprint = argc - 1;
280 ROL::Ptr<std::ostream> outStream;
283 outStream = ROL::makePtrFromRef(std::cout);
285 outStream = ROL::makePtrFromRef(bhs);
299 ROL::Ptr<std::vector<RealT>> x_ptr = ROL::makePtr<std::vector<RealT>>(
dim, 0.0);
300 ROL::Ptr<std::vector<RealT>> g_ptr = ROL::makePtr<std::vector<RealT>>(
dim, 0.0);
301 ROL::Ptr<std::vector<RealT>> l_ptr = ROL::makePtr<std::vector<RealT>>(
dim, 0.0);
303 for (
int i=0; i<dim/2; i++) {
304 (*x_ptr)[2*i] = -1.2;
305 (*x_ptr)[2*i+1] = 1.0;
308 (*l_ptr)[2*i] = ROL::ROL_NINF<RealT>();
309 (*l_ptr)[2*i+1] = -1.5;
311 ROL::Ptr<OptStdVector<RealT>> x = ROL::makePtr<OptStdVector<RealT>>(x_ptr);
312 ROL::Ptr<OptDualStdVector<RealT>> g = ROL::makePtr<OptDualStdVector<RealT>>(g_ptr);
313 ROL::Ptr<OptStdVector<RealT>> l = ROL::makePtr<OptStdVector<RealT>>(l_ptr);
316 ROL::Ptr<ROL::Bounds<RealT>> bnd = ROL::makePtr<ROL::Bounds<RealT>>(*l);
319 ROL::Ptr<std::vector<RealT> > aa_ptr = ROL::makePtr<std::vector<RealT>>(1, 1.0);
321 ROL::Ptr<std::vector<RealT> > bb_ptr = ROL::makePtr<std::vector<RealT>>(1, 2.0);
323 ROL::Ptr<std::vector<RealT> > cc_ptr = ROL::makePtr<std::vector<RealT>>(1, 3.0);
325 std::vector<RealT> std_vec_err = av.
checkVector(bv,cv,
true,*outStream);
328 ROL::Ptr<ROL::Problem<RealT>> problem
329 = ROL::makePtr<ROL::Problem<RealT>>(obj,x,g);
330 problem->addBoundConstraint(bnd);
331 problem->finalize(
false,
true,*outStream);
334 ROL::ParameterList parlist;
335 std::string stepname =
"Trust Region";
336 parlist.sublist(
"Step").set(
"Type",stepname);
338 parlist.sublist(
"Step").sublist(stepname).set(
"Subproblem Solver",
"Truncated CG");
339 parlist.sublist(
"Step").sublist(stepname).set(
"Subproblem Model",
"Lin-More");
340 parlist.sublist(
"General").set(
"Output Level",1);
341 parlist.sublist(
"General").sublist(
"Krylov").set(
"Relative Tolerance",1e-2);
342 parlist.sublist(
"General").sublist(
"Krylov").set(
"Iteration Limit",10);
343 parlist.sublist(
"General").sublist(
"Krylov").set(
"Absolute Tolerance",1e-4);
344 parlist.sublist(
"General").sublist(
"Secant").set(
"Type",
"Limited-Memory BFGS");
345 parlist.sublist(
"General").sublist(
"Secant").set(
"Use as Hessian",
false);
346 parlist.sublist(
"Status Test").set(
"Gradient Tolerance",1.e-12);
347 parlist.sublist(
"Status Test").set(
"Step Tolerance",1.e-14);
348 parlist.sublist(
"Status Test").set(
"Iteration Limit",100);
349 ROL::Ptr<ROL::Solver<RealT>> solver
350 = ROL::makePtr<ROL::Solver<RealT>>(problem,parlist);
353 solver->solve(*outStream);
356 ROL::Ptr<std::vector<RealT> > xtrue_ptr = ROL::makePtr<std::vector<RealT>>(
dim, 1.0);
360 x->axpy(-1.0, xtrue);
361 RealT abserr = x->norm();
363 *outStream << std::scientific <<
"\n Absolute solution error: " << abserr;
364 *outStream << std::scientific <<
"\n Relative solution error: " << relerr;
365 if ( relerr > sqrt(ROL::ROL_EPSILON<RealT>()) ) {
368 ROL::Ptr<std::vector<RealT> > vec_err_ptr = ROL::makePtr<std::vector<RealT>>(std_vec_err);
370 *outStream << std::scientific <<
"\n Linear algebra error: " << vec_err.norm() << std::endl;
371 if ( vec_err.norm() > 1e2*ROL::ROL_EPSILON<RealT>() ) {
375 catch (std::logic_error& err) {
376 *outStream << err.what() <<
"\n";
381 std::cout <<
"End Result: TEST FAILED\n";
383 std::cout <<
"End Result: TEST PASSED\n";
typename PV< Real >::size_type size_type
void scale(const Real alpha)
Compute where .
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< std::vector< Element > > getVector()
OptDualStdVector(const ROL::Ptr< std::vector< Element > > &std_vec)
virtual std::vector< Real > checkVector(const Vector< Real > &x, const Vector< Real > &y, const bool printToStream=true, std::ostream &outStream=std::cout) const
Verify vector-space methods.
Contains definitions for Rosenbrock's function.
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...
void applyUnary(const ROL::Elementwise::UnaryFunction< Real > &f)
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.
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
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.
void applyUnary(const ROL::Elementwise::UnaryFunction< Real > &f)
Real reduce(const ROL::Elementwise::ReductionOp< Real > &r) const
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
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.
Real reduce(const ROL::Elementwise::ReductionOp< Real > &r) const
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 .
void applyBinary(const ROL::Elementwise::BinaryFunction< Real > &f, const ROL::Vector< Real > &x)
Real norm() const
Returns where .
int dimension() const
Return dimension of the vector space.
ROL::Ptr< std::vector< Element > > std_vec_
int main(int argc, char *argv[])
void applyBinary(const ROL::Elementwise::BinaryFunction< Real > &f, const ROL::Vector< Real > &x)
ROL::Ptr< const std::vector< Element > > getVector() const
Real dot(const ROL::Vector< Real > &x) const
Compute where .
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.