55 #include "Teuchos_GlobalMPISession.hpp"
65 template <
class Real,
class Element=Real>
68 template <
class Real,
class Element=Real>
75 template <
class Real,
class Element>
85 mutable ROL::Ptr<OptDualStdVector<Real> >
dual_vec_;
96 (*std_vec_)[i] += (*xvalptr)[i];
103 (*std_vec_)[i] *= alpha;
113 val += (*std_vec_)[i]*(*xvalptr)[i];
120 val = std::sqrt(
dot(*
this) );
124 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
125 return ROL::makePtr<OptStdVector>( ROL::makePtr<std::vector<Element>>(
std_vec_->size()) );
128 ROL::Ptr<const std::vector<Element> >
getVector()
const {
136 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
138 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(
std_vec_->size(),0.0);
141 ROL::Ptr<V> e = ROL::makePtr<OptStdVector>( e_ptr );
149 dual_vec_ = ROL::makePtr<OptDualStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
157 template <
class Real,
class Element>
179 (*std_vec_)[i] += (*xvalptr)[i];
186 (*std_vec_)[i] *= alpha;
196 val += (*std_vec_)[i]*(*xvalptr)[i];
203 val = std::sqrt(
dot(*
this) );
207 ROL::Ptr<ROL::Vector<Real> >
clone()
const {
208 return ROL::makePtr<OptDualStdVector>( ROL::makePtr<std::vector<Element>>(
std_vec_->size()) );
211 ROL::Ptr<const std::vector<Element> >
getVector()
const {
219 ROL::Ptr<ROL::Vector<Real> >
basis(
const int i )
const {
221 ROL::Ptr<vector> e_ptr = ROL::makePtr<vector>(
std_vec_->size(),0.0);
224 ROL::Ptr<V> e = ROL::makePtr<OptDualStdVector>( e_ptr );
231 dual_vec_ = ROL::makePtr<OptStdVector<Real>>( ROL::makePtr<std::vector<Element>>(*std_vec_) );
245 int main(
int argc,
char *argv[]) {
247 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
250 int iprint = argc - 1;
251 ROL::Ptr<std::ostream> outStream;
254 outStream = ROL::makePtrFromRef(std::cout);
256 outStream = ROL::makePtrFromRef(bhs);
268 ROL::ParameterList parlist;
269 std::string stepname =
"Trust Region";
270 parlist.sublist(
"Step").sublist(stepname).set(
"Subproblem Solver",
"Truncated CG");
271 parlist.sublist(
"General").sublist(
"Krylov").set(
"Iteration Limit",10);
272 parlist.sublist(
"General").sublist(
"Krylov").set(
"Relative Tolerance",1e-2);
273 parlist.sublist(
"General").sublist(
"Krylov").set(
"Absolute Tolerance",1e-4);
274 parlist.sublist(
"General").sublist(
"Secant").set(
"Use as Hessian",
true);
275 parlist.sublist(
"Status Test").set(
"Gradient Tolerance",1.e-12);
276 parlist.sublist(
"Status Test").set(
"Step Tolerance",1.e-14);
277 parlist.sublist(
"Status Test").set(
"Iteration Limit",100);
278 ROL::Ptr<ROL::Step<RealT>>
279 step = ROL::makePtr<ROL::TrustRegionStep<RealT>>(parlist);
280 ROL::Ptr<ROL::StatusTest<RealT>>
281 status = ROL::makePtr<ROL::StatusTest<RealT>>(parlist);
285 ROL::Ptr<std::vector<RealT> > x_ptr = ROL::makePtr<std::vector<RealT>>(
dim, 0.0);
286 ROL::Ptr<std::vector<RealT> > g_ptr = ROL::makePtr<std::vector<RealT>>(
dim, 0.0);
288 for (
int i=0; i<dim/2; i++) {
289 (*x_ptr)[2*i] = -1.2;
290 (*x_ptr)[2*i+1] = 1.0;
298 ROL::Ptr<std::vector<RealT> > aa_ptr = ROL::makePtr<std::vector<RealT>>(1, 1.0);
300 ROL::Ptr<std::vector<RealT> > bb_ptr = ROL::makePtr<std::vector<RealT>>(1, 2.0);
302 ROL::Ptr<std::vector<RealT> > cc_ptr = ROL::makePtr<std::vector<RealT>>(1, 3.0);
304 std::vector<RealT> std_vec_err = av.
checkVector(bv,cv,
true,*outStream);
307 std::vector<std::string> output = algo.run(x,g, obj,
true, *outStream);
310 ROL::Ptr<std::vector<RealT> > xtrue_ptr = ROL::makePtr<std::vector<RealT>>(
dim, 1.0);
317 *outStream << std::scientific <<
"\n Absolute solution error: " << abserr;
318 *outStream << std::scientific <<
"\n Relative solution error: " << relerr;
319 if ( relerr > sqrt(ROL::ROL_EPSILON<RealT>()) ) {
322 ROL::Ptr<std::vector<RealT> > vec_err_ptr = ROL::makePtr<std::vector<RealT>>(std_vec_err);
324 *outStream << std::scientific <<
"\n Linear algebra error: " << vec_err.norm() << std::endl;
325 if ( vec_err.norm() > 1e2*ROL::ROL_EPSILON<RealT>() ) {
329 catch (std::logic_error& err) {
330 *outStream << err.what() <<
"\n";
335 std::cout <<
"End Result: TEST FAILED\n";
337 std::cout <<
"End Result: TEST PASSED\n";
typename PV< Real >::size_type size_type
Rosenbrock's function.
void scale(const Real alpha)
Compute where .
std::vector< Element > vector
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
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...
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.
Provides an interface to run optimization algorithms.
ROL::Ptr< OptDualStdVector< Real > > dual_vec_
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.