ROL
burgers-control/example_03.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 
15 #include "example_03.hpp"
16 
17 typedef double RealT;
18 
19 int main(int argc, char *argv[]) {
20 
21  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
22 
23  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
24  int iprint = argc - 1;
25  ROL::Ptr<std::ostream> outStream;
26  ROL::nullstream bhs; // outputs nothing
27  if (iprint > 0)
28  outStream = ROL::makePtrFromRef(std::cout);
29  else
30  outStream = ROL::makePtrFromRef(bhs);
31 
32  int errorFlag = 0;
33 
34  // *** Example body.
35 
36  try {
37  // Initialize full objective function.
38  int nx = 80; // Set spatial discretization.
39  int nt = 80; // Set temporal discretization.
40  RealT T = 1.0; // Set end time.
41  RealT alpha = 5e-2; // Set penalty parameter.
42  RealT nu = 1e-2; // Set viscosity parameter.
43  Objective_BurgersControl<RealT> obj(alpha,nx,nt,T);
44  // Initialize equality constraints
45  Constraint_BurgersControl<RealT> con(nx, nt, T, nu);
46  // Initialize iteration vectors.
47  ROL::Ptr<std::vector<RealT> > z_ptr = ROL::makePtr<std::vector<RealT>>((nx+2)*(nt+1), 1.0);
48  ROL::Ptr<std::vector<RealT> > gz_ptr = ROL::makePtr<std::vector<RealT>>((nx+2)*(nt+1), 1.0);
49  ROL::Ptr<std::vector<RealT> > yz_ptr = ROL::makePtr<std::vector<RealT>>((nx+2)*(nt+1), 1.0);
50  for (int i=0; i<(nx+2)*(nt+1); i++) {
51  (*z_ptr)[i] = (RealT)rand()/(RealT)RAND_MAX;
52  (*yz_ptr)[i] = (RealT)rand()/(RealT)RAND_MAX;
53  }
54  ROL::StdVector<RealT> z(z_ptr);
55  ROL::StdVector<RealT> gz(gz_ptr);
56  ROL::StdVector<RealT> yz(yz_ptr);
57  ROL::Ptr<ROL::Vector<RealT> > zp = ROL::makePtrFromRef(z);
58  ROL::Ptr<ROL::Vector<RealT> > gzp = ROL::makePtrFromRef(gz);
59  ROL::Ptr<ROL::Vector<RealT> > yzp = ROL::makePtrFromRef(yz);
60 
61  ROL::Ptr<std::vector<RealT> > u_ptr = ROL::makePtr<std::vector<RealT>>(nx*nt, 1.0);
62  ROL::Ptr<std::vector<RealT> > gu_ptr = ROL::makePtr<std::vector<RealT>>(nx*nt, 1.0);
63  ROL::Ptr<std::vector<RealT> > yu_ptr = ROL::makePtr<std::vector<RealT>>(nx*nt, 1.0);
64  for (int i=0; i<nx*nt; i++) {
65  (*u_ptr)[i] = (RealT)rand()/(RealT)RAND_MAX;
66  (*yu_ptr)[i] = (RealT)rand()/(RealT)RAND_MAX;
67  }
68  ROL::StdVector<RealT> u(u_ptr);
69  ROL::StdVector<RealT> gu(gu_ptr);
70  ROL::StdVector<RealT> yu(yu_ptr);
71  ROL::Ptr<ROL::Vector<RealT> > up = ROL::makePtrFromRef(u);
72  ROL::Ptr<ROL::Vector<RealT> > gup = ROL::makePtrFromRef(gu);
73  ROL::Ptr<ROL::Vector<RealT> > yup = ROL::makePtrFromRef(yu);
74 
75  ROL::Ptr<std::vector<RealT> > c_ptr = ROL::makePtr<std::vector<RealT>>(nx*nt, 1.0);
76  ROL::Ptr<std::vector<RealT> > l_ptr = ROL::makePtr<std::vector<RealT>>(nx*nt, 1.0);
77  ROL::StdVector<RealT> c(c_ptr);
78  ROL::StdVector<RealT> l(l_ptr);
79 
81  ROL::Vector_SimOpt<RealT> g(gup,gzp);
82  ROL::Vector_SimOpt<RealT> y(yup,yzp);
83  // Check derivatives.
84  obj.checkGradient(x,x,y,true,*outStream);
85  obj.checkHessVec(x,x,y,true,*outStream);
86  con.checkApplyJacobian(x,y,c,true,*outStream);
87  //con.checkApplyAdjointJacobian(x,yu,c,x,true,*outStream);
88  con.checkApplyAdjointHessian(x,yu,y,x,true,*outStream);
89  // Check consistency of Jacobians and adjoint Jacobians.
90  con.checkAdjointConsistencyJacobian_1(c,yu,u,z,true,*outStream);
91  con.checkAdjointConsistencyJacobian_2(c,yz,u,z,true,*outStream);
92  // Check consistency of solves.
93  con.checkSolve(u,z,c,true,*outStream);
94  con.checkInverseJacobian_1(c,yu,u,z,true,*outStream);
95  con.checkInverseAdjointJacobian_1(yu,c,u,z,true,*outStream);
96 
97  // Initialize reduced objective function.
98  ROL::Ptr<std::vector<RealT> > p_ptr = ROL::makePtr<std::vector<RealT>>(nx*nt, 1.0);
99  ROL::StdVector<RealT> p(p_ptr);
100  ROL::Ptr<ROL::Vector<RealT> > pp = ROL::makePtrFromRef(p);
101  ROL::Ptr<ROL::Objective_SimOpt<RealT> > pobj = ROL::makePtrFromRef(obj);
102  ROL::Ptr<ROL::Constraint_SimOpt<RealT> > pcon = ROL::makePtrFromRef(con);
103  ROL::Reduced_Objective_SimOpt<RealT> robj(pobj,pcon,up,zp,pp);
104  // Check derivatives.
105  robj.checkGradient(z,z,yz,true,*outStream);
106  robj.checkHessVec(z,z,yz,true,*outStream);
107  // Get input parameter list.
108  std::string filename = "input.xml";
109  auto parlist = ROL::getParametersFromXmlFile( filename );
110  parlist->sublist("Status Test").set("Gradient Tolerance",1.e-10);
111  parlist->sublist("Status Test").set("Constraint Tolerance",1.e-10);
112  parlist->sublist("Status Test").set("Step Tolerance",1.e-16);
113  parlist->sublist("Status Test").set("Iteration Limit",100);
114  // Build Algorithm pointer.
115  ROL::Ptr<ROL::Algorithm<RealT>> algo;
116  ROL::Ptr<ROL::Step<RealT>> step;
117  ROL::Ptr<ROL::StatusTest<RealT>> status;
118 
119  // Solve using trust regions.
120  step = ROL::makePtr<ROL::TrustRegionStep<RealT>>(*parlist);
121  status = ROL::makePtr<ROL::StatusTest<RealT>>(*parlist);
122  algo = ROL::makePtr<ROL::Algorithm<RealT>>(step,status,false);
123  z.zero();
124  std::clock_t timer_tr = std::clock();
125  algo->run(z,robj,true,*outStream);
126  *outStream << "Trust-Region Newton required " << (std::clock()-timer_tr)/(RealT)CLOCKS_PER_SEC
127  << " seconds.\n";
128  ROL::Ptr<ROL::Vector<RealT> > zTR = z.clone();
129  zTR->set(z);
130 
131  // Solve using a composite step method.
132  step = ROL::makePtr<ROL::CompositeStep<RealT>>(*parlist);
133  status = ROL::makePtr<ROL::ConstraintStatusTest<RealT>>(*parlist);
134  algo = ROL::makePtr<ROL::Algorithm<RealT>>(step,status,false);
135  x.zero();
136  ROL::Elementwise::Fill<RealT> setFunc(0.25);
137  x.applyUnary(setFunc);
138  std::clock_t timer_cs = std::clock();
139  algo->run(x,g,l,c,obj,con,true,*outStream);
140  *outStream << "Composite Step required " << (std::clock()-timer_cs)/(RealT)CLOCKS_PER_SEC
141  << " seconds.\n";
142 
143  // Compute error between solutions
144  ROL::Ptr<ROL::Vector<RealT> > err = z.clone();
145  err->set(*zTR); err->axpy(-1.,z);
146  errorFlag += (err->norm() > 1.e-4) ? 1 : 0;
147  if (errorFlag) {
148  *outStream << "\n\nControl error = " << err->norm() << "\n";
149  }
150 
151 // std::ofstream control;
152 // control.open("control.txt");
153 // for (int t = 0; t < nt+1; t++) {
154 // for (int n = 0; n < nx+2; n++) {
155 // control << (RealT)t/(RealT)nt << " "
156 // << (RealT)n/((RealT)(nx+1)) << " "
157 // << (*z_ptr)[t*(nx+2)+n] << "\n";
158 // }
159 // }
160 // control.close();
161 //
162 // std::ofstream state;
163 // state.open("state.txt");
164 // for (int t = 0; t < nt; t++) {
165 // for (int n = 0; n < nx; n++) {
166 // state << (RealT)(t+1)/(RealT)nt << " "
167 // << (RealT)(n+1)/((RealT)(nx+1)) << " "
168 // << (*u_ptr)[t*nx+n] << "\n";
169 // }
170 // }
171 // state.close();
172  }
173  catch (std::logic_error& err) {
174  *outStream << err.what() << "\n";
175  errorFlag = -1000;
176  }; // end try
177 
178  if (errorFlag != 0)
179  std::cout << "End Result: TEST FAILED\n";
180  else
181  std::cout << "End Result: TEST PASSED\n";
182 
183  return 0;
184 
185 }
186 
virtual Real checkSolve(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &c, const bool printToStream=true, std::ostream &outStream=std::cout)
Defines the linear algebra or vector space interface for simulation-based optimization.
virtual Real checkAdjointConsistencyJacobian_2(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the primary interface.
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:133
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.
basic_nullstream< char, std::char_traits< char >> nullstream
Definition: ROL_Stream.hpp:36
void applyUnary(const Elementwise::UnaryFunction< Real > &f)
virtual Real checkInverseJacobian_1(const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
virtual Ptr< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
virtual std::vector< std::vector< Real > > checkApplyAdjointHessian(const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &step, const bool printToScreen=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the application of the adjoint of constraint Hessian. ...
virtual std::vector< std::vector< Real > > checkApplyJacobian(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the constraint Jacobian application.
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.
virtual Real checkAdjointConsistencyJacobian_1(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the primary interface.
virtual Real checkInverseAdjointJacobian_1(const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)