ROL
poisson-control/example_02.cpp
Go to the documentation of this file.
1 
2 // Burgers includes
3 #include "example_02.hpp"
4 // ROL includes
5 #include "ROL_Algorithm.hpp"
6 #include "ROL_StdVector.hpp"
7 #include "ROL_StdTeuchosBatchManager.hpp"
11 #include "ROL_Vector_SimOpt.hpp"
12 #include "ROL_Bounds.hpp"
13 #include "ROL_ParameterList.hpp"
14 
15 // Teuchos includes
16 #include "Teuchos_Time.hpp"
17 #include "ROL_Stream.hpp"
18 #include "Teuchos_GlobalMPISession.hpp"
19 #include "Teuchos_Comm.hpp"
20 #include "Teuchos_DefaultComm.hpp"
21 #include "Teuchos_CommHelpers.hpp"
22 
23 int main( int argc, char *argv[] ) {
24 
25  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
26 
27  auto comm = ROL::toPtr( Teuchos::DefaultComm<int>::getComm() );
28 
29  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
30  int iprint = argc - 1;
31  ROL::Ptr<std::ostream> outStream;
32  ROL::nullstream bhs; // outputs nothing
33  if (iprint > 0)
34  outStream = ROL::makePtrFromRef(std::cout);
35  else
36  outStream = ROL::makePtrFromRef(bhs);
37 
38  int errorFlag = 0;
39 
40  // *** Example body.
41 
42  try {
43 
44  /***************************************************************************/
45  /***************** GRAB INPUTS *********************************************/
46  /***************************************************************************/
47  // Get finite element parameter list
48  std::string filename = "example_02.xml";
49  auto parlist = ROL::getParametersFromXmlFile( filename );
50 
51  if ( parlist->get("Display Option",0) && (comm->getRank() > 0) ) {
52  parlist->set("Display Option",0);
53  }
54  // Get ROL parameter list
55  filename = "input.xml";
56  auto ROL_parlist = ROL::getParametersFromXmlFile( filename );
57 
58  /***************************************************************************/
59  /***************** INITIALIZE SAMPLERS *************************************/
60  /***************************************************************************/
61  int dim = 2;
62  bool useSA = parlist->get("Use Stochastic Approximation",false);
63  int nSamp = 1;
64  if ( !useSA ) {
65  nSamp = parlist->get("Number of Monte Carlo Samples",1000);
66  }
67  std::vector<double> tmp(2); tmp[0] = -1.0; tmp[1] = 1.0;
68  std::vector<std::vector<double> > bounds(dim,tmp);
69  ROL::Ptr<ROL::BatchManager<double> > bman
70  = ROL::makePtr<ROL::StdTeuchosBatchManager<double,int>>(comm);
71  ROL::Ptr<ROL::SampleGenerator<double> > sampler
72  = ROL::makePtr<ROL::MonteCarloGenerator<double>>(nSamp,bounds,bman,useSA);
73 
74  /***************************************************************************/
75  /***************** INITIALIZE CONTROL VECTOR *******************************/
76  /***************************************************************************/
77  int nx = parlist->get("Number of Elements", 128);
78  ROL::Ptr<std::vector<double> > z_ptr = ROL::makePtr<std::vector<double>>(nx+1, 0.0);
79  ROL::Ptr<ROL::Vector<double> > z = ROL::makePtr<ROL::StdVector<double>>(z_ptr);
80  ROL::Ptr<std::vector<double> > u_ptr = ROL::makePtr<std::vector<double>>(nx-1, 0.0);
81  ROL::Ptr<ROL::Vector<double> > u = ROL::makePtr<ROL::StdVector<double>>(u_ptr);
83  ROL::Ptr<std::vector<double> > p_ptr = ROL::makePtr<std::vector<double>>(nx-1, 0.0);
84  ROL::Ptr<ROL::Vector<double> > p = ROL::makePtr<ROL::StdVector<double>>(p_ptr);
85  ROL::Ptr<std::vector<double> > U_ptr = ROL::makePtr<std::vector<double>>(nx+1, 35.0);
86  ROL::Ptr<ROL::Vector<double> > U = ROL::makePtr<ROL::StdVector<double>>(U_ptr);
87  ROL::Ptr<std::vector<double> > L_ptr = ROL::makePtr<std::vector<double>>(nx+1, -5.0);
88  ROL::Ptr<ROL::Vector<double> > L = ROL::makePtr<ROL::StdVector<double>>(L_ptr);
89  ROL::Bounds<double> bnd(L,U);
90 
91  /***************************************************************************/
92  /***************** INITIALIZE OBJECTIVE FUNCTION ***************************/
93  /***************************************************************************/
94  double alpha = parlist->get("Penalty Parameter", 1.e-4);
95  ROL::Ptr<FEM<double> > fem = ROL::makePtr<FEM<double>>(nx);
96  ROL::Ptr<ROL::Objective_SimOpt<double> > pObj
97  = ROL::makePtr<DiffusionObjective<double>>(fem, alpha);
98  ROL::Ptr<ROL::Constraint_SimOpt<double> > pCon
99  = ROL::makePtr<DiffusionConstraint<double>>(fem);
100  ROL::Ptr<ROL::Objective<double> > robj
101  = ROL::makePtr<ROL::Reduced_Objective_SimOpt<double>>(pObj,pCon,u,z,p);
102  ROL::RiskNeutralObjective<double> obj(robj,sampler);
103 
104  /***************************************************************************/
105  /***************** RUN DERIVATIVE CHECK ************************************/
106  /***************************************************************************/
107  if (parlist->get("Run Derivative Check",false)) {
108  // Direction to test finite differences
109  ROL::Ptr<std::vector<double> > dz_ptr = ROL::makePtr<std::vector<double>>(nx+1, 0.0);
110  ROL::Ptr<ROL::Vector<double> > dz = ROL::makePtr<ROL::StdVector<double>>(dz_ptr);
111  ROL::Ptr<std::vector<double> > du_ptr = ROL::makePtr<std::vector<double>>(nx-1, 0.0);
112  ROL::Ptr<ROL::Vector<double> > du = ROL::makePtr<ROL::StdVector<double>>(du_ptr);
114  // Set to random vectors
115  srand(12345);
116  for (int i=0; i<nx+1; i++) {
117  (*dz_ptr)[i] = 2.0*(double)rand()/(double)RAND_MAX - 1.0;
118  (*z_ptr)[i] = 2.0*(double)rand()/(double)RAND_MAX - 1.0;
119  }
120  for (int i=0; i<nx-1; i++) {
121  (*du_ptr)[i] = 2.0*(double)rand()/(double)RAND_MAX - 1.0;
122  (*u_ptr)[i] = 2.0*(double)rand()/(double)RAND_MAX - 1.0;
123  }
124  // Run derivative checks
125  std::vector<double> param(dim,0.0);
126  robj->setParameter(param);
127  if ( comm->getRank() == 0 ) {
128  std::cout << "\nRUN DERIVATIVE CHECK FOR PARAMETRIZED OBJECTIVE FUNCTION SIMOPT\n";
129  }
130  pObj->checkGradient(x,d,(comm->getRank()==0));
131  pObj->checkHessVec(x,d,(comm->getRank()==0));
132  if ( comm->getRank() == 0 ) {
133  std::cout << "\nRUN DERIVATIVE CHECK FOR PARAMETRIZED EQUALITY CONSTRAINT SIMOPT\n";
134  }
135  pCon->checkApplyJacobian(x,d,*p,(comm->getRank()==0));
136  pCon->checkApplyAdjointJacobian(x,*du,*p,x,(comm->getRank()==0));
137  pCon->checkApplyAdjointHessian(x,*du,d,x,(comm->getRank()==0));
138  if ( comm->getRank() == 0 ) {
139  std::cout << "\nRUN DERIVATIVE CHECK FOR PARAMETRIZED OBJECTIVE FUNCTION\n";
140  }
141  robj->checkGradient(*z,*dz,(comm->getRank()==0));
142  robj->checkHessVec(*z,*dz,(comm->getRank()==0));
143  // Run derivative checks
144  if ( comm->getRank() == 0 ) {
145  std::cout << "\nRUN DERIVATIVE CHECK FOR RISK-NEUTRAL OBJECTIVE FUNCTION\n";
146  }
147  obj.checkGradient(*z,*dz,(comm->getRank()==0));
148  obj.checkHessVec(*z,*dz,(comm->getRank()==0));
149  }
150 
151  /***************************************************************************/
152  /***************** INITIALIZE ROL ALGORITHM ********************************/
153  /***************************************************************************/
154  ROL::Ptr<ROL::Algorithm<double> > algo;
155  if ( useSA ) {
156  ROL_parlist->sublist("General").set("Recompute Objective Function",false);
157  ROL_parlist->sublist("Step").sublist("Line Search").set("Initial Step Size",0.1/alpha);
158  ROL_parlist->sublist("Step").sublist("Line Search").set("User Defined Initial Step Size",true);
159  ROL_parlist->sublist("Step").sublist("Line Search").sublist("Line-Search Method").set("Type","Iteration Scaling");
160  ROL_parlist->sublist("Step").sublist("Line Search").sublist("Descent Method").set("Type","Steepest Descent");
161  ROL_parlist->sublist("Step").sublist("Line Search").sublist("Curvature Condition").set("Type","Null Curvature Condition");
162  algo = ROL::makePtr<ROL::Algorithm<double>>("Line Search",*ROL_parlist,false);
163  }
164  else {
165  algo = ROL::makePtr<ROL::Algorithm<double>>("Trust Region",*ROL_parlist,false);
166  }
167 
168  /***************************************************************************/
169  /***************** PERFORM OPTIMIZATION ************************************/
170  /***************************************************************************/
171  Teuchos::Time timer("Optimization Time",true);
172  z->zero();
173  algo->run(*z,obj,bnd,(comm->getRank()==0));
174  double optTime = timer.stop();
175 
176  /***************************************************************************/
177  /***************** PRINT RESULTS *******************************************/
178  /***************************************************************************/
179  int my_number_samples = sampler->numMySamples(), number_samples = 0;
180  Teuchos::reduceAll<int,int>(*comm,Teuchos::REDUCE_SUM,1,&my_number_samples,&number_samples);
181  int my_number_solves = ROL::dynamicPtrCast<DiffusionConstraint<double> >(pCon)->getNumSolves(), number_solves = 0;
182  Teuchos::reduceAll<int,int>(*comm,Teuchos::REDUCE_SUM,1,&my_number_solves,&number_solves);
183  if (comm->getRank() == 0) {
184  std::cout << "Number of Samples = " << number_samples << "\n";
185  std::cout << "Number of Solves = " << number_solves << "\n";
186  std::cout << "Optimization Time = " << optTime << "\n\n";
187  }
188 
189  if ( comm->getRank() == 0 ) {
190  std::ofstream file;
191  if (useSA) {
192  file.open("control_SA.txt");
193  }
194  else {
195  file.open("control_SAA.txt");
196  }
197  std::vector<double> xmesh(fem->nz(),0.0);
198  fem->build_mesh(xmesh);
199  for (int i = 0; i < fem->nz(); i++ ) {
200  file << std::setprecision(std::numeric_limits<double>::digits10) << std::scientific << xmesh[i] << " "
201  << std::setprecision(std::numeric_limits<double>::digits10) << std::scientific << (*z_ptr)[i]
202  << "\n";
203  }
204  file.close();
205  }
206  }
207  catch (std::logic_error err) {
208  *outStream << err.what() << "\n";
209  errorFlag = -1000;
210  }; // end try
211 
212  if (errorFlag != 0)
213  std::cout << "End Result: TEST FAILED\n";
214  else
215  std::cout << "End Result: TEST PASSED\n";
216 
217  return 0;
218 }
219 
220 
221 
222 
Defines the linear algebra or vector space interface for simulation-based optimization.
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
Provides the elementwise interface to apply upper and lower bound constraints.
Definition: ROL_Bounds.hpp:59
basic_nullstream< char, char_traits< char >> nullstream
Definition: ROL_Stream.hpp:72
int main(int argc, char *argv[])