53 #include "Teuchos_GlobalMPISession.hpp"
63 template <
class Real,
class Element=Real>
66 template <
class Real,
class Element=Real>
73 template <
class Real,
class Element>
83 mutable ROL::Ptr<OptDualStdVector<Real> >
dual_vec_;
94 (*std_vec_)[i] += (*xvalptr)[i];
98 void scale(
const Real alpha ) {
101 (*std_vec_)[i] *= alpha;
111 val += (*std_vec_)[i]*(*xvalptr)[i];
118 val = std::sqrt(
dot(*
this) );
122 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
123 return ROL::makePtr<OptStdVector>( ROL::makePtr<std::vector<Element>>(
std_vec_->size()) );
126 ROL::Ptr<const std::vector<Element> >
getVector()
const {
134 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
136 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(
std_vec_->size(),0.0);
139 ROL::Ptr<V> e = ROL::makePtr<OptStdVector>( e_ptr );
147 dual_vec_ = ROL::makePtr<OptDualStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
157 val += (*std_vec_)[i]*(*xvalptr)[i];
166 template <
class Real,
class Element>
188 (*std_vec_)[i] += (*xvalptr)[i];
195 (*std_vec_)[i] *= alpha;
205 val += (*std_vec_)[i]*(*xvalptr)[i];
212 val = std::sqrt(
dot(*
this) );
216 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
217 return ROL::makePtr<OptDualStdVector>( ROL::makePtr<std::vector<Element>>(
std_vec_->size()) );
220 ROL::Ptr<const std::vector<Element> >
getVector()
const {
228 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
230 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(
std_vec_->size(),0.0);
233 ROL::Ptr<V> e = ROL::makePtr<OptDualStdVector>( e_ptr );
240 dual_vec_ = ROL::makePtr<OptStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
250 val += (*std_vec_)[i]*(*xvalptr)[i];
265 int main(
int argc,
char *argv[]) {
267 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
270 int iprint = argc - 1;
271 ROL::Ptr<std::ostream> outStream;
274 outStream = ROL::makePtrFromRef(std::cout);
276 outStream = ROL::makePtrFromRef(bhs);
290 ROL::Ptr<std::vector<RealT>> x_ptr = ROL::makePtr<std::vector<RealT>>(
dim, 0.0);
291 ROL::Ptr<std::vector<RealT>> g_ptr = ROL::makePtr<std::vector<RealT>>(
dim, 0.0);
293 for (
int i=0; i<dim/2; i++) {
294 (*x_ptr)[2*i] = -1.2;
295 (*x_ptr)[2*i+1] = 1.0;
299 ROL::Ptr<OptStdVector<RealT>> x = ROL::makePtr<OptStdVector<RealT>>(x_ptr);
300 ROL::Ptr<OptDualStdVector<RealT>> g = ROL::makePtr<OptDualStdVector<RealT>>(g_ptr);
303 ROL::Ptr<std::vector<RealT> > aa_ptr = ROL::makePtr<std::vector<RealT>>(1, 1.0);
305 ROL::Ptr<std::vector<RealT> > bb_ptr = ROL::makePtr<std::vector<RealT>>(1, 2.0);
307 ROL::Ptr<std::vector<RealT> > cc_ptr = ROL::makePtr<std::vector<RealT>>(1, 3.0);
309 std::vector<RealT> std_vec_err = av.
checkVector(bv,cv,
true,*outStream);
312 ROL::Ptr<ROL::Problem<RealT>> problem
313 = ROL::makePtr<ROL::Problem<RealT>>(obj,x,g);
314 problem->finalize(
false,
true,*outStream);
317 ROL::ParameterList parlist;
318 std::string stepname =
"Trust Region";
319 parlist.sublist(
"Step").set(
"Type",stepname);
321 parlist.sublist(
"Step").sublist(stepname).set(
"Subproblem Solver",
"Truncated CG");
322 parlist.sublist(
"General").set(
"Output Level",1);
323 parlist.sublist(
"General").sublist(
"Krylov").set(
"Relative Tolerance",1e-2);
324 parlist.sublist(
"General").sublist(
"Krylov").set(
"Iteration Limit",10);
325 parlist.sublist(
"General").sublist(
"Krylov").set(
"Absolute Tolerance",1e-4);
326 parlist.sublist(
"General").sublist(
"Secant").set(
"Type",
"Limited-Memory BFGS");
327 parlist.sublist(
"General").sublist(
"Secant").set(
"Use as Hessian",
true);
328 parlist.sublist(
"Status Test").set(
"Gradient Tolerance",1.e-12);
329 parlist.sublist(
"Status Test").set(
"Step Tolerance",1.e-14);
330 parlist.sublist(
"Status Test").set(
"Iteration Limit",100);
331 ROL::Ptr<ROL::Solver<RealT>> solver
332 = ROL::makePtr<ROL::Solver<RealT>>(problem,parlist);
335 solver->solve(*outStream);
338 ROL::Ptr<std::vector<RealT> > xtrue_ptr = ROL::makePtr<std::vector<RealT>>(
dim, 1.0);
342 x->axpy(-1.0, xtrue);
343 RealT abserr = x->norm();
345 *outStream << std::scientific <<
"\n Absolute solution error: " << abserr;
346 *outStream << std::scientific <<
"\n Relative solution error: " << relerr;
347 if ( relerr > sqrt(ROL::ROL_EPSILON<RealT>()) ) {
350 ROL::Ptr<std::vector<RealT> > vec_err_ptr = ROL::makePtr<std::vector<RealT>>(std_vec_err);
352 *outStream << std::scientific <<
"\n Linear algebra error: " << vec_err.norm() << std::endl;
353 if ( vec_err.norm() > 1e2*ROL::ROL_EPSILON<RealT>() ) {
357 catch (std::logic_error& err) {
358 *outStream << err.what() <<
"\n";
363 std::cout <<
"End Result: TEST FAILED\n";
365 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...
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()
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.
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.
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 .
Real norm() const
Returns where .
int dimension() const
Return dimension of the vector space.
ROL::Ptr< std::vector< Element > > std_vec_
basic_nullstream< char, char_traits< char >> nullstream
int main(int argc, char *argv[])
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.