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