ROL
sacado/example_03.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 
55 #include <iostream>
56 
57 
58 #include "ROL_Sacado_Objective_SimOpt.hpp"
59 #include "ROL_Sacado_EqualityConstraint_SimOpt.hpp"
60 #include "ROL_Vector_SimOpt.hpp"
61 
62 #include "ROL_LineSearchStep.hpp"
63 #include "ROL_Algorithm.hpp"
65 #include "ROL_CompositeStepSQP.hpp"
66 
67 #include "Teuchos_oblackholestream.hpp"
68 #include "Teuchos_GlobalMPISession.hpp"
69 #include "Teuchos_XMLParameterListHelpers.hpp"
70 
71 #include "example_03.hpp"
72 
73 using namespace ROL;
74 
75 typedef double RealT;
76 
77 int main(int argc, char **argv)
78 {
79 
80 
81  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
82 
83  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
84  int iprint = argc - 1;
85  Teuchos::RCP<std::ostream> outStream;
86  Teuchos::oblackholestream bhs; // outputs nothing
87  if (iprint > 0)
88  outStream = Teuchos::rcp(&std::cout, false);
89  else
90  outStream = Teuchos::rcp(&bhs, false);
91 
92  int errorFlag = 0;
93 
94  // *** Example body.
95 
96  try {
97 
98  Teuchos::ParameterList parlist;
99  parlist.set("Nominal SQP Optimality Solver Tolerance", 1.e-2);
100  ROL::CompositeStepSQP<RealT> step(parlist);
101 
102  int ni = 12; // Number of interpolation points
103  int nq = 20; // Number of quadrature points
104 
105  RealT xl = 0.0; // Left end of domain
106  RealT xr = 5.0; // Right end of domain
107 
108  Teuchos::RCP<NodalBasis<RealT> > basisp = Teuchos::rcp(new NodalBasis<RealT> (ni,nq) );
109 
110  // simulation vector (state u) and optimization vector (control z)
111  Teuchos::RCP<std::vector<RealT> > u_rcp = Teuchos::rcp(new std::vector<RealT>(ni,0));
112  Teuchos::RCP<std::vector<RealT> > z_rcp = Teuchos::rcp(new std::vector<RealT>(ni,0));
113 
114  // change of simulation vector (state u) and optimization vector (control z)
115  Teuchos::RCP<std::vector<RealT> > yu_rcp = Teuchos::rcp(new std::vector<RealT>(ni,0));
116  Teuchos::RCP<std::vector<RealT> > yz_rcp = Teuchos::rcp(new std::vector<RealT>(ni,0));
117 
118  // Gradients
119  Teuchos::RCP<std::vector<RealT> > gu_rcp = Teuchos::rcp(new std::vector<RealT>(ni,0));
120  Teuchos::RCP<std::vector<RealT> > gz_rcp = Teuchos::rcp(new std::vector<RealT>(ni,0));
121 
122  // Constraint and Lagrange multiplier value
123  Teuchos::RCP<std::vector<RealT> > c_rcp = Teuchos::rcp(new std::vector<RealT>(ni,1.0));
124  Teuchos::RCP<std::vector<RealT> > l_rcp = Teuchos::rcp(new std::vector<RealT>(ni,1.0));
125 
126 
127 
128  for(int i=0; i<ni; ++i) {
129  (*u_rcp)[i] = (RealT)rand()/(RealT)RAND_MAX;
130  (*z_rcp)[i] = (RealT)rand()/(RealT)RAND_MAX;
131  (*yu_rcp)[i] = (RealT)rand()/(RealT)RAND_MAX;
132  (*yz_rcp)[i] = (RealT)rand()/(RealT)RAND_MAX;
133 
134  }
135 
136 
137  StdVector<RealT> u(u_rcp);
138  StdVector<RealT> z(z_rcp);
139  StdVector<RealT> yu(yu_rcp);
140  StdVector<RealT> yz(yz_rcp);
141  StdVector<RealT> gu(gu_rcp);
142  StdVector<RealT> gz(gz_rcp);
143  StdVector<RealT> c(c_rcp);
144  StdVector<RealT> l(l_rcp);
145 
146  Teuchos::RCP<Vector<RealT> > zp = Teuchos::rcp(&z,false);
147  Teuchos::RCP<Vector<RealT> > up = Teuchos::rcp(&u,false);
148  Teuchos::RCP<Vector<RealT> > yzp = Teuchos::rcp(&yz,false);
149  Teuchos::RCP<Vector<RealT> > yup = Teuchos::rcp(&yu,false);
150  Teuchos::RCP<Vector<RealT> > gzp = Teuchos::rcp(&gz,false);
151  Teuchos::RCP<Vector<RealT> > gup = Teuchos::rcp(&gu,false);
152 
153  // Set tracking target
154  Teuchos::RCP<std::vector<RealT> > u_targ_rcp = Teuchos::rcp(new std::vector<RealT>(ni,0));
155 
156  for(int i=0;i<ni;++i) {
157  (*u_targ_rcp)[i] = tanh(2.0*((*basisp->xip_)[i]-2.5));
158  }
159 
160  // Regularization parameter
161  RealT gamma = 1e-4;
162 
163  StdVector<RealT> u_targ(u_targ_rcp);
164 
165  BoundaryValueProblem<RealT> bvp(xl,xr,basisp);
166  QuadraticTracking<RealT> quadtrack(gamma,u_targ,basisp);
167 
168  // Instantiate constraint and objective objects
170  Sacado_Objective_SimOpt<RealT,QuadraticTracking> obj(quadtrack);
171 
172 
173  Vector_SimOpt<RealT> x(up,zp);
174  Vector_SimOpt<RealT> g(gup,gzp);
175  Vector_SimOpt<RealT> y(yup,yzp);
176 
177  // Check derivatives
178  obj.checkGradient(x,x,y,true,*outStream);
179  obj.checkHessVec(x,x,y,true,*outStream);
180  constr.checkApplyJacobian(x,y,c,true,*outStream);
181  constr.checkApplyAdjointJacobian(x,yu,c,x,true,*outStream);
182  constr.checkApplyAdjointHessian(x,yu,y,x,true,*outStream);
183  try {
184  constr.checkInverseJacobian_1(c,yu,u,z,true,*outStream);
185  } catch (const std::logic_error &e) {
186  *outStream << e.what();
187  }
188  try {
189  constr.checkInverseAdjointJacobian_1(c,yu,u,z,true,*outStream);
190  } catch (const std::logic_error &e) {
191  *outStream << e.what();
192  }
193  try {
194  constr.checkSolve(u,z,c,true,*outStream);
195  } catch (const std::logic_error &e) {
196  *outStream << e.what();
197  }
198 
199  // Define Status Test
200  RealT gtol = 1e-12; // norm of gradient tolerance
201  RealT ctol = 1e-12; // norm of constraint tolerance
202  RealT stol = 1e-18; // norm of step tolerance
203  int maxit = 1000; // maximum number of iterations
204  ROL::StatusTestSQP<RealT> status(gtol, ctol, stol, maxit);
205 
206  // Define Algorithm
207  ROL::DefaultAlgorithm<RealT> algo(step, status, false);
208  algo.run(x,g,l,c,obj,constr,true,*outStream);
209 
210  }
211  catch (std::logic_error err) {
212  *outStream << err.what() << "\n";
213  errorFlag = -1000;
214  }; // end try
215 
216  if (errorFlag != 0)
217  std::cout << "End Result: TEST FAILED\n";
218  else
219  std::cout << "End Result: TEST PASSED\n";
220 
221  return 0;
222 
223 
224 }
Defines the linear algebra or vector space interface for simulation-based optimization.
Inherit and add method for applying the inverse partial constraint Jacobian and its adjoint...
Provides the std::vector implementation of the ROL::Vector interface.
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.
Implements the computation of optimization steps with composite-step trust-region SQP methods...
Objective to minimize .
Compute the constraint vector -u'' + exp(u) - z with natural boundary conditions. ...
double RealT
int main(int argc, char **argv)
double RealT