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