ROL
dual-spaces/rosenbrock-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 
48 #define USE_HESSVEC 1
49 
50 #include "ROL_Rosenbrock.hpp"
51 #include "ROL_LineSearchStep.hpp"
52 #include "ROL_TrustRegionStep.hpp"
53 #include "ROL_Algorithm.hpp"
54 #include "Teuchos_oblackholestream.hpp"
55 #include "Teuchos_GlobalMPISession.hpp"
56 #include <iostream>
57 
58 typedef double RealT;
59 
60 
61 /*** Declare two vector spaces. ***/
62 
63 // Forward declarations:
64 
65 template <class Real, class Element=Real>
66 class OptStdVector; // Optimization space.
67 
68 template <class Real, class Element=Real>
69 class OptDualStdVector; // Dual optimization space.
70 
71 
72 // Vector space definitions:
73 
74 // Optimization space.
75 template <class Real, class Element>
76 class OptStdVector : public ROL::Vector<Real> {
77 
78 private:
79 Teuchos::RCP<std::vector<Element> > std_vec_;
80 mutable Teuchos::RCP<OptDualStdVector<Real> > dual_vec_;
81 
82 public:
83 
84 OptStdVector(const Teuchos::RCP<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(Teuchos::null) {}
85 
86 void plus( const ROL::Vector<Real> &x ) {
87  OptStdVector &ex = Teuchos::dyn_cast<OptStdVector>(const_cast <ROL::Vector<Real> &>(x));
88  Teuchos::RCP<const std::vector<Element> > xvalptr = ex.getVector();
89  unsigned dimension = std_vec_->size();
90  for (unsigned i=0; i<dimension; i++) {
91  (*std_vec_)[i] += (*xvalptr)[i];
92  }
93 }
94 
95 void scale( const Real alpha ) {
96  unsigned dimension = std_vec_->size();
97  for (unsigned i=0; i<dimension; i++) {
98  (*std_vec_)[i] *= alpha;
99  }
100 }
101 
102 Real dot( const ROL::Vector<Real> &x ) const {
103  Real val = 0;
104  OptStdVector<Real, Element> & ex = Teuchos::dyn_cast<OptStdVector<Real, Element> >(const_cast <ROL::Vector<Real> &>(x));
105  Teuchos::RCP<const std::vector<Element> > xvalptr = ex.getVector();
106  unsigned dimension = std_vec_->size();
107  for (unsigned i=0; i<dimension; i++) {
108  val += (*std_vec_)[i]*(*xvalptr)[i];
109  }
110  return val;
111 }
112 
113 Real norm() const {
114  Real val = 0;
115  val = std::sqrt( dot(*this) );
116  return val;
117 }
118 
119 Teuchos::RCP<ROL::Vector<Real> > clone() const {
120  return Teuchos::rcp( new OptStdVector( Teuchos::rcp( new std::vector<Element>(std_vec_->size()) ) ) );
121 }
122 
123 Teuchos::RCP<const std::vector<Element> > getVector() const {
124  return std_vec_;
125 }
126 
127 Teuchos::RCP<ROL::Vector<Real> > basis( const int i ) const {
128  Teuchos::RCP<OptStdVector> e = Teuchos::rcp( new OptStdVector( Teuchos::rcp(new std::vector<Element>(std_vec_->size(), 0.0)) ) );
129  (const_cast <std::vector<Element> &> (*e->getVector()))[i]= 1.0;
130  return e;
131 }
132 
133 int dimension() const {return std_vec_->size();}
134 
135 const ROL::Vector<Real> & dual() const {
136  dual_vec_ = Teuchos::rcp( new OptDualStdVector<Real>( Teuchos::rcp( new std::vector<Element>(*std_vec_) ) ) );
137  return *dual_vec_;
138 }
139 
140 }; // class OptStdVector
141 
142 
143 // Dual optimization space.
144 template <class Real, class Element>
145 class OptDualStdVector : public ROL::Vector<Real> {
146 
147 private:
148 Teuchos::RCP<std::vector<Element> > std_vec_;
149 mutable Teuchos::RCP<OptStdVector<Real> > dual_vec_;
150 
151 public:
152 
153 OptDualStdVector(const Teuchos::RCP<std::vector<Element> > & std_vec) : std_vec_(std_vec), dual_vec_(Teuchos::null) {}
154 
155 void plus( const ROL::Vector<Real> &x ) {
156  OptDualStdVector &ex = Teuchos::dyn_cast<OptDualStdVector>(const_cast <ROL::Vector<Real> &>(x));
157  Teuchos::RCP<const std::vector<Element> > xvalptr = ex.getVector();
158  unsigned dimension = std_vec_->size();
159  for (unsigned i=0; i<dimension; i++) {
160  (*std_vec_)[i] += (*xvalptr)[i];
161  }
162 }
163 
164 void scale( const Real alpha ) {
165  unsigned dimension = std_vec_->size();
166  for (unsigned i=0; i<dimension; i++) {
167  (*std_vec_)[i] *= alpha;
168  }
169 }
170 
171 Real dot( const ROL::Vector<Real> &x ) const {
172  Real val = 0;
173  OptDualStdVector<Real, Element> & ex = Teuchos::dyn_cast<OptDualStdVector<Real, Element> >(const_cast <ROL::Vector<Real> &>(x));
174  Teuchos::RCP<const std::vector<Element> > xvalptr = ex.getVector();
175  unsigned dimension = std_vec_->size();
176  for (unsigned i=0; i<dimension; i++) {
177  val += (*std_vec_)[i]*(*xvalptr)[i];
178  }
179  return val;
180 }
181 
182 Real norm() const {
183  Real val = 0;
184  val = std::sqrt( dot(*this) );
185  return val;
186 }
187 
188 Teuchos::RCP<ROL::Vector<Real> > clone() const {
189  return Teuchos::rcp( new OptDualStdVector( Teuchos::rcp( new std::vector<Element>(std_vec_->size()) ) ) );
190 }
191 
192 Teuchos::RCP<const std::vector<Element> > getVector() const {
193  return std_vec_;
194 }
195 
196 Teuchos::RCP<ROL::Vector<Real> > basis( const int i ) const {
197  Teuchos::RCP<OptDualStdVector> e = Teuchos::rcp( new OptDualStdVector( Teuchos::rcp(new std::vector<Element>(std_vec_->size(), 0.0)) ) );
198  (const_cast <std::vector<Element> &> (*e->getVector()))[i]= 1.0;
199  return e;
200 }
201 
202 int dimension() const {return std_vec_->size();}
203 
204 const ROL::Vector<Real> & dual() const {
205  dual_vec_ = Teuchos::rcp( new OptStdVector<Real>( Teuchos::rcp( new std::vector<Element>(*std_vec_) ) ) );
206  return *dual_vec_;
207 }
208 
209 }; // class OptDualStdVector
210 
211 
212 /*** End of declaration of two vector spaces. ***/
213 
214 
215 
216 
217 
218 
219 int main(int argc, char *argv[]) {
220 
221  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
222 
223  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
224  int iprint = argc - 1;
225  Teuchos::RCP<std::ostream> outStream;
226  Teuchos::oblackholestream bhs; // outputs nothing
227  if (iprint > 0)
228  outStream = Teuchos::rcp(&std::cout, false);
229  else
230  outStream = Teuchos::rcp(&bhs, false);
231 
232  int errorFlag = 0;
233 
234  // *** Example body.
235 
236  try {
237 
239  int dim = 100; // Set problem dimension. Must be even.
240 
241  Teuchos::ParameterList parlist;
242  // Enumerations
243  parlist.set("Descent Type", "Quasi-Newton Method");
244  parlist.set("Secant Type", "Limited-memory BFGS");
245  parlist.set("Linesearch Type", "Cubic Interpolation");
246  parlist.set("Linesearch Curvature Condition", "Wolfe");
247  // Linesearch Parameters
248  parlist.set("Maximum Number of Function Evaluations", 20);
249  parlist.set("Sufficient Decrease Parameter", 1.e-4);
250  parlist.set("Curvature Conditions Parameter", 0.9);
251  parlist.set("Backtracking Rate", 0.5);
252  parlist.set("Initial Linesearch Parameter", 1.0);
253  parlist.set("User Defined Linesearch Parameter", false);
254  // Krylov Parameters
255  parlist.set("Absolute Krylov Tolerance", 1.e-4);
256  parlist.set("Relative Krylov Tolerance", 1.e-2);
257  parlist.set("Maximum Number of Krylov Iterations", 10);
258  // Trust Region Parameters
259  parlist.set("Trust-Region Subproblem Solver Type","Truncated CG");
260  parlist.set("Use Secant Hessian-Times-A-Vector",true);
261  // Define Step
262  //ROL::LineSearchStep<RealT> step(parlist);
263  ROL::TrustRegionStep<RealT> step(parlist);
264 
265 
266  // Define Status Test
267  RealT gtol = 1e-12; // norm of gradient tolerance
268  RealT stol = 1e-14; // norm of step tolerance
269  int maxit = 100; // maximum number of iterations
270  ROL::StatusTest<RealT> status(gtol, stol, maxit);
271 
272  // Define Algorithm
273  ROL::DefaultAlgorithm<RealT> algo(step,status,false);
274 
275  // Iteration Vector
276  Teuchos::RCP<std::vector<RealT> > x_rcp = Teuchos::rcp( new std::vector<RealT> (dim, 0.0) );
277  Teuchos::RCP<std::vector<RealT> > g_rcp = Teuchos::rcp( new std::vector<RealT> (dim, 0.0) );
278  // Set Initial Guess
279  for (int i=0; i<dim/2; i++) {
280  (*x_rcp)[2*i] = -1.2;
281  (*x_rcp)[2*i+1] = 1.0;
282  (*g_rcp)[2*i] = 0;
283  (*g_rcp)[2*i+1] = 0;
284  }
285 
286 
287  OptStdVector<RealT> x(x_rcp); // Iteration Vector
288  OptDualStdVector<RealT> g(g_rcp); // zeroed gradient vector in dual space
289 
290  Teuchos::RCP<std::vector<RealT> > aa_rcp = Teuchos::rcp( new std::vector<RealT> (1, 1.0) );
291  OptDualStdVector<RealT> av(aa_rcp);
292  Teuchos::RCP<std::vector<RealT> > bb_rcp = Teuchos::rcp( new std::vector<RealT> (1, 2.0) );
293  OptDualStdVector<RealT> bv(bb_rcp);
294  Teuchos::RCP<std::vector<RealT> > cc_rcp = Teuchos::rcp( new std::vector<RealT> (1, 3.0) );
295  OptDualStdVector<RealT> cv(cc_rcp);
296  av.checkVector(bv,cv);
297 
298 
299  // Run Algorithm
300  std::vector<std::string> output = algo.run(x,g, obj, false);
301  for ( unsigned i = 0; i < output.size(); i++ ) {
302  std::cout << output[i];
303  }
304 
305  // Get True Solution
306  Teuchos::RCP<std::vector<RealT> > xtrue_rcp = Teuchos::rcp( new std::vector<RealT> (dim, 1.0) );
307  OptStdVector<RealT> xtrue(xtrue_rcp);
308 
309  // Compute Error
310  x.axpy(-1.0, xtrue);
311  RealT abserr = x.norm();
312  RealT relerr = abserr/xtrue.norm();
313  *outStream << std::scientific << "\n Absolute Error: " << abserr;
314  *outStream << std::scientific << "\n Relative Error: " << relerr << "\n";
315  if ( relerr > sqrt(ROL::ROL_EPSILON) ) {
316  errorFlag += 1;
317  }
318  }
319  catch (std::logic_error err) {
320  *outStream << err.what() << "\n";
321  errorFlag = -1000;
322  }; // end try
323 
324  if (errorFlag != 0)
325  std::cout << "End Result: TEST FAILED\n";
326  else
327  std::cout << "End Result: TEST PASSED\n";
328 
329  return 0;
330 
331 }
332 
Teuchos::RCP< std::vector< Element > > std_vec_
Teuchos::RCP< OptDualStdVector< Real > > dual_vec_
Rosenbrock's function.
void scale(const Real alpha)
Compute where .
OptStdVector(const Teuchos::RCP< std::vector< Element > > &std_vec)
Teuchos::RCP< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:141
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
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.
Definition: ROL_Vector.hpp:242
int main(int argc, char *argv[])
Teuchos::RCP< OptStdVector< Real > > dual_vec_
Teuchos::RCP< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Contains definitions for Rosenbrock's function.
Teuchos::RCP< const std::vector< Element > > getVector() const
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
OptDualStdVector(const Teuchos::RCP< std::vector< Element > > &std_vec)
Teuchos::RCP< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
void scale(const Real alpha)
Compute where .
virtual std::vector< std::string > run(Vector< Real > &x, Objective< Real > &obj, bool print=false, std::ostream &outStream=std::cout)
Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface.
Real dot(const ROL::Vector< Real > &x) const
Compute where .
int dimension() const
Return dimension of the vector space.
Real norm() const
Returns where .
Provides an interface to check status of optimization algorithms.
Teuchos::RCP< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
int dimension() const
Return dimension of the vector space.
Real dot(const ROL::Vector< Real > &x) const
Compute where .
double RealT
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:194
Teuchos::RCP< const std::vector< Element > > getVector() const
Teuchos::RCP< std::vector< Element > > std_vec_
void plus(const ROL::Vector< Real > &x)
Compute , where .
void plus(const ROL::Vector< Real > &x)
Compute , where .
Real norm() const
Returns where .
Provides the interface to compute optimization steps with trust regions.
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:115