ROL
function/test_05.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Rapid Optimization Library (ROL) Package
4 //
5 // Copyright 2014 NTESS and the ROL contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
17 #include "ROL_StdVector.hpp"
19 #include "ROL_Algorithm.hpp"
20 #include "ROL_TrustRegionStep.hpp"
21 #include "ROL_Stream.hpp"
22 #include "Teuchos_GlobalMPISession.hpp"
23 
24 #include <iostream>
25 
26 typedef double RealT;
27 
28 
29 int main(int argc, char *argv[]) {
30 
31  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
32 
33  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
34  int iprint = argc - 1;
35  ROL::Ptr<std::ostream> outStream;
36  ROL::nullstream bhs; // outputs nothing
37  if (iprint > 0)
38  outStream = ROL::makePtrFromRef(std::cout);
39  else
40  outStream = ROL::makePtrFromRef(bhs);
41 
42  int errorFlag = 0;
43 
44  // *** Example body.
45 
46  try {
47 
48  ROL::Ptr<ROL::Objective<RealT>> obj;
49  ROL::Ptr<ROL::Constraint<RealT>> constr;
50  ROL::Ptr<ROL::Vector<RealT>> x;
51  ROL::Ptr<ROL::Vector<RealT>> sol;
52 
53  // Retrieve objective, constraint, iteration vector, solution vector.
55  obj = SEC.getObjective();
56  constr = SEC.getEqualityConstraint();
57  x = SEC.getInitialGuess();
58  sol = SEC.getSolution();
59 
60  // Inititalize vectors
61  int dim = 5;
62  int nc = 3;
63  RealT left = -1e0, right = 1e0;
64  ROL::Ptr<std::vector<RealT>> xtest_ptr = ROL::makePtr<std::vector<RealT>>(dim, 0.0);
65  ROL::Ptr<std::vector<RealT>> g_ptr = ROL::makePtr<std::vector<RealT>>(dim, 0.0);
66  ROL::Ptr<std::vector<RealT>> d_ptr = ROL::makePtr<std::vector<RealT>>(dim, 0.0);
67  ROL::Ptr<std::vector<RealT>> v_ptr = ROL::makePtr<std::vector<RealT>>(dim, 0.0);
68  ROL::Ptr<std::vector<RealT>> vc_ptr = ROL::makePtr<std::vector<RealT>>(nc, 0.0);
69  ROL::Ptr<std::vector<RealT>> vl_ptr = ROL::makePtr<std::vector<RealT>>(nc, 0.0);
70  ROL::StdVector<RealT> xtest(xtest_ptr);
71  ROL::StdVector<RealT> g(g_ptr);
72  ROL::StdVector<RealT> d(d_ptr);
73  ROL::StdVector<RealT> v(v_ptr);
74  ROL::StdVector<RealT> vc(vc_ptr);
75  ROL::StdVector<RealT> vl(vl_ptr);
76  // set xtest, d, v
77  for (int i=0; i<dim; i++) {
78  (*xtest_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
79  (*d_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
80  (*v_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
81  }
82  // set vc, vl
83  for (int i=0; i<nc; i++) {
84  (*vc_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
85  (*vl_ptr)[i] = ( (RealT)rand() / (RealT)RAND_MAX ) * (right - left) + left;
86  }
87 
88  xtest.set(*x);
89 
90  // Initialize nonlinear least squares objectives
91  ROL::NonlinearLeastSquaresObjective<RealT> nlls(constr,*x,vc,false);
92  ROL::NonlinearLeastSquaresObjective<RealT> gnnlls(constr,*x,vc,true);
93 
94  // Check derivatives
95  constr->checkApplyJacobian(xtest, v, vc, true, *outStream); *outStream << "\n";
96  constr->checkApplyAdjointJacobian(xtest, vl, vc, xtest, true, *outStream); *outStream << "\n";
97  constr->checkApplyAdjointHessian(xtest, vl, d, xtest, true, *outStream); *outStream << "\n";
98  nlls.checkGradient(xtest, d, true, *outStream); *outStream << "\n";
99  nlls.checkHessVec(xtest, v, true, *outStream); *outStream << "\n";
100  nlls.checkHessSym(xtest, d, v, true, *outStream); *outStream << "\n";
101 
102  // Define algorithm.
103  ROL::ParameterList parlist;
104  parlist.sublist("Step").sublist("Trust Region").set("Subproblem Solver","Truncated CG");
105  parlist.sublist("Status Test").set("Gradient Tolerance",1.e-10);
106  parlist.sublist("Status Test").set("Constraint Tolerance",1.e-10);
107  parlist.sublist("Status Test").set("Step Tolerance",1.e-18);
108  parlist.sublist("Status Test").set("Iteration Limit",100);
109  ROL::Ptr<ROL::StatusTest<RealT>> status = ROL::makePtr<ROL::StatusTest<RealT>>(parlist);
110  ROL::Ptr<ROL::Step<RealT>> step = ROL::makePtr<ROL::TrustRegionStep<RealT>>(parlist);
111  ROL::Algorithm<RealT> algo(step,status,false);
112 
113  // Run Algorithm
114  *outStream << "\nSOLVE USING FULL HESSIAN\n";
115  x->set(xtest);
116  algo.run(*x, nlls, true, *outStream);
117  algo.reset();
118  *outStream << "\nSOLVE USING GAUSS-NEWTON HESSIAN\n";
119  x->set(xtest);
120  algo.run(*x, gnnlls, true, *outStream);
121  }
122  catch (std::logic_error& err) {
123  *outStream << err.what() << "\n";
124  errorFlag = -1000;
125  }; // end try
126 
127  if (errorFlag != 0)
128  std::cout << "End Result: TEST FAILED\n";
129  else
130  std::cout << "End Result: TEST PASSED\n";
131 
132  return 0;
133 
134 }
135 
Ptr< Constraint< Real > > getEqualityConstraint(void) const
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
virtual std::vector< std::vector< Real > > checkGradient(const Vector< Real > &x, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference gradient check.
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
Provides an interface to run optimization algorithms.
Ptr< Objective< Real > > getObjective(void) const
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.
void set(const Vector< Real > &x)
Set where .
Ptr< Vector< Real > > getSolution(const int i=0) const
basic_nullstream< char, char_traits< char >> nullstream
Definition: ROL_Stream.hpp:38
int main(int argc, char *argv[])
virtual std::vector< std::vector< Real > > checkHessVec(const Vector< Real > &x, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference Hessian-applied-to-vector check.
Provides the interface to evaluate nonlinear least squares objective functions.
virtual std::vector< Real > checkHessSym(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout)
Hessian symmetry check.
constexpr auto dim
Ptr< Vector< Real > > getInitialGuess(void) const