22 #include "Teuchos_GlobalMPISession.hpp"
26 int main(
int argc,
char *argv[]) {
28 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
32 int iprint = argc - 1;
33 ROL::Ptr<std::ostream> outStream;
36 outStream = ROL::makePtrFromRef(std::cout);
38 outStream = ROL::makePtrFromRef(bhs);
43 RealT tol = std::sqrt(ROL::ROL_EPSILON<RealT>());
45 ROL::Ptr<ROL::Vector<RealT>> sol, mul, x, lam, l, u, c;
46 ROL::Ptr<ROL::Objective<RealT>> obj;
47 ROL::Ptr<ROL::Constraint<RealT>> con;
48 ROL::Ptr<ROL::BoundConstraint<RealT>> bnd;
49 ROL::Ptr<ROL::PolyhedralProjection<RealT>> proj;
50 ROL::ParameterList list;
51 list.sublist(
"General").set(
"Output Level",2);
52 std::vector<RealT> data;
54 *outStream << std::endl <<
"Hock and Schittkowski Problem #41" << std::endl << std::endl;
62 lam = mul->clone(); lam->set(*mul);
63 x = sol->clone(); x->set(*sol);
64 l = sol->clone(); l->zero();
65 u = sol->clone(); u->setScalar(static_cast<RealT>(1));
66 c = mul->dual().clone();
68 list.sublist(
"General").sublist(
"Polyhedral Projection").set(
"Type",
"Dai-Fletcher");
69 proj = ROL::PolyhedralProjectionFactory<RealT>(*sol,sol->dual(),bnd,con,*lam,*c,list);
70 proj->project(*x,*outStream);
72 con->value(*c,*x,tol);
75 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(sol)->getVector();
76 *outStream <<
" Initial: x1 = " << data[0] <<
" x2 = " << data[1]
77 <<
" x3 = " << data[2] << std::endl;
78 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(x)->getVector();
79 *outStream <<
" Result: x1 = " << data[0] <<
" x2 = " << data[1]
80 <<
" x3 = " << data[2] << std::endl;
81 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(lam)->getVector();
82 *outStream <<
" Multiplier: l1 = " << data[0] << std::endl;
84 *outStream << std::endl;
85 *outStream <<
" is equality feasible = " << (cnorm<=tol) << std::endl
86 <<
" are bounds feasible = " << bnd->isFeasible(*x) << std::endl;
88 errorFlag += !bnd->isFeasible(*x);
89 errorFlag += (cnorm > tol);
91 *outStream << std::endl <<
"Hock and Schittkowski Problem #41" << std::endl << std::endl;
99 lam = mul->clone(); lam->set(*mul);
100 x = sol->clone(); x->set(*sol);
101 l = sol->clone(); l->zero();
102 u = sol->clone(); u->setScalar(static_cast<RealT>(1));
103 c = mul->dual().clone();
105 list.sublist(
"General").sublist(
"Polyhedral Projection").set(
"Type",
"Ridders");
106 proj = ROL::PolyhedralProjectionFactory<RealT>(*sol,sol->dual(),bnd,con,*lam,*c,list);
107 proj->project(*x,*outStream);
109 con->value(*c,*x,tol);
112 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(sol)->getVector();
113 *outStream <<
" Initial: x1 = " << data[0] <<
" x2 = " << data[1]
114 <<
" x3 = " << data[2] << std::endl;
115 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(x)->getVector();
116 *outStream <<
" Result: x1 = " << data[0] <<
" x2 = " << data[1]
117 <<
" x3 = " << data[2] << std::endl;
118 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(lam)->getVector();
119 *outStream <<
" Multiplier: l1 = " << data[0] << std::endl;
121 *outStream << std::endl;
122 *outStream <<
" is equality feasible = " << (cnorm<=tol) << std::endl
123 <<
" are bounds feasible = " << bnd->isFeasible(*x) << std::endl;
125 errorFlag += !bnd->isFeasible(*x);
126 errorFlag += (cnorm > tol);
128 *outStream << std::endl <<
"Hock and Schittkowski Problem #41" << std::endl << std::endl;
136 lam = mul->clone(); lam->set(*mul);
137 x = sol->clone(); x->set(*sol);
138 l = sol->clone(); l->zero();
139 u = sol->clone(); u->setScalar(static_cast<RealT>(1));
140 c = mul->dual().clone();
142 list.sublist(
"General").sublist(
"Polyhedral Projection").set(
"Type",
"Brents");
143 proj = ROL::PolyhedralProjectionFactory<RealT>(*sol,sol->dual(),bnd,con,*lam,*c,list);
144 proj->project(*x,*outStream);
146 con->value(*c,*x,tol);
149 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(sol)->getVector();
150 *outStream <<
" Initial: x1 = " << data[0] <<
" x2 = " << data[1]
151 <<
" x3 = " << data[2] << std::endl;
152 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(x)->getVector();
153 *outStream <<
" Result: x1 = " << data[0] <<
" x2 = " << data[1]
154 <<
" x3 = " << data[2] << std::endl;
155 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(lam)->getVector();
156 *outStream <<
" Multiplier: l1 = " << data[0] << std::endl;
158 *outStream << std::endl;
159 *outStream <<
" is equality feasible = " << (cnorm<=tol) << std::endl
160 <<
" are bounds feasible = " << bnd->isFeasible(*x) << std::endl;
162 errorFlag += !bnd->isFeasible(*x);
163 errorFlag += (cnorm > tol);
165 *outStream << std::endl <<
"Hock and Schittkowski Problem #53" << std::endl << std::endl;
173 lam = mul->clone(); lam->set(*mul);
174 x = sol->clone(); x->set(*sol);
175 l = sol->clone(); l->zero();
176 u = sol->clone(); u->setScalar(static_cast<RealT>(1));
177 c = mul->dual().clone();
179 list.sublist(
"General").sublist(
"Polyhedral Projection").set(
"Type",
"Dykstra");
180 proj = ROL::PolyhedralProjectionFactory<RealT>(*sol,sol->dual(),bnd,con,*lam,*c,list);
181 proj->project(*x,*outStream);
183 con->value(*c,*x,tol);
186 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(sol)->getVector();
187 *outStream <<
" Initial: x1 = " << data[0] <<
" x2 = " << data[1]
188 <<
" x3 = " << data[2] <<
" x4 = " << data[3]
189 <<
" x5 = " << data[4] << std::endl;
190 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(x)->getVector();
191 *outStream <<
" Result: x1 = " << data[0] <<
" x2 = " << data[1]
192 <<
" x3 = " << data[2] <<
" x4 = " << data[3]
193 <<
" x5 = " << data[4] << std::endl;
194 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(lam)->getVector();
195 *outStream <<
" Multiplier: l1 = " << data[0] <<
" l2 = " << data[1]
196 <<
" l3 = " << data[2] << std::endl;
198 *outStream << std::endl;
199 *outStream <<
" is equality feasible = " << (cnorm<=tol) << std::endl
200 <<
" are bounds feasible = " << bnd->isFeasible(*x) << std::endl;
202 errorFlag += !bnd->isFeasible(*x);
203 errorFlag += (cnorm > tol);
205 *outStream << std::endl <<
"Hock and Schittkowski Problem #53" << std::endl << std::endl;
213 lam = mul->clone(); lam->set(*mul);
214 x = sol->clone(); x->set(*sol);
215 l = sol->clone(); l->zero();
216 u = sol->clone(); u->setScalar(static_cast<RealT>(1));
217 c = mul->dual().clone();
219 list.sublist(
"General").sublist(
"Polyhedral Projection").set(
"Type",
"Douglas-Rachford");
220 proj = ROL::PolyhedralProjectionFactory<RealT>(*sol,sol->dual(),bnd,con,*lam,*c,list);
221 proj->project(*x,*outStream);
223 con->value(*c,*x,tol);
226 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(sol)->getVector();
227 *outStream <<
" Initial: x1 = " << data[0] <<
" x2 = " << data[1]
228 <<
" x3 = " << data[2] <<
" x4 = " << data[3]
229 <<
" x5 = " << data[4] << std::endl;
230 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(x)->getVector();
231 *outStream <<
" Result: x1 = " << data[0] <<
" x2 = " << data[1]
232 <<
" x3 = " << data[2] <<
" x4 = " << data[3]
233 <<
" x5 = " << data[4] << std::endl;
234 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(lam)->getVector();
235 *outStream <<
" Multiplier: l1 = " << data[0] <<
" l2 = " << data[1]
236 <<
" l3 = " << data[2] << std::endl;
238 *outStream << std::endl;
239 *outStream <<
" is equality feasible = " << (cnorm<=tol) << std::endl
240 <<
" are bounds feasible = " << bnd->isFeasible(*x) << std::endl;
242 errorFlag += !bnd->isFeasible(*x);
243 errorFlag += (cnorm > tol);
245 *outStream << std::endl <<
"Hock and Schittkowski Problem #55" << std::endl << std::endl;
259 lam = mul->clone(); lam->set(*mul);
260 x = sol->clone(); x->set(*sol);
261 l = sol->clone(); l->zero();
262 u = sol->clone(); u->setScalar(static_cast<RealT>(1));
263 c = mul->dual().clone();
265 list.sublist(
"General").sublist(
"Polyhedral Projection").set(
"Type",
"Semismooth Newton");
266 proj = ROL::PolyhedralProjectionFactory<RealT>(*sol,sol->dual(),bnd,con,*lam,*c,list);
267 proj->project(*x,*outStream);
269 con->value(*c,*x,tol);
272 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(sol)->getVector();
273 *outStream <<
" Initial: x1 = " << data[0] <<
" x2 = " << data[1]
274 <<
" x3 = " << data[2] <<
" x4 = " << data[3]
275 <<
" x5 = " << data[4] <<
" x6 = " << data[5] << std::endl;
276 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(x)->getVector();
277 *outStream <<
" Result: x1 = " << data[0] <<
" x2 = " << data[1]
278 <<
" x3 = " << data[2] <<
" x4 = " << data[3]
279 <<
" x5 = " << data[4] <<
" x6 = " << data[5] << std::endl;
280 data = *ROL::staticPtrCast<ROL::StdVector<RealT>>(lam)->getVector();
281 *outStream <<
" Multiplier: l1 = " << data[0] <<
" l2 = " << data[1]
282 <<
" l3 = " << data[2] <<
" l4 = " << data[3]
283 <<
" l5 = " << data[4] <<
" l6 = " << data[5] << std::endl;
285 *outStream << std::endl;
286 *outStream <<
" is equality feasible = " << (cnorm<=tol) << std::endl
287 <<
" are bounds feasible = " << bnd->isFeasible(*x) << std::endl;
288 *outStream << std::endl;
290 errorFlag += !bnd->isFeasible(*x);
291 errorFlag += (cnorm > tol);
294 catch (std::logic_error& err) {
295 *outStream << err.what() <<
"\n";
300 std::cout <<
"End Result: TEST FAILED\n";
302 std::cout <<
"End Result: TEST PASSED\n";
Ptr< Vector< Real > > getInitialGuess(void) const
Contains definitions for W. Hock and K. Schittkowski 55th test function.
Ptr< BoundConstraint< Real > > getBoundConstraint(void) const
Ptr< Constraint< Real > > getEqualityConstraint(void) const
Contains definitions for W. Hock and K. Schittkowski 41th test function.
Contains definitions for W. Hock and K. Schittkowski 53th test function.
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
Ptr< Constraint< Real > > getEqualityConstraint(void) const
basic_nullstream< char, std::char_traits< char >> nullstream
Ptr< BoundConstraint< Real > > getBoundConstraint(void) const
Ptr< Objective< Real > > getObjective(void) const
Ptr< Vector< Real > > getEqualityMultiplier(void) const
Ptr< Vector< Real > > getEqualityMultiplier(void) const
Ptr< Vector< Real > > getInitialGuess(void) const
Ptr< Objective< Real > > getObjective(void) const
int main(int argc, char *argv[])
Ptr< Objective< Real > > getObjective(void) const
Ptr< BoundConstraint< Real > > getBoundConstraint(void) const
Ptr< Vector< Real > > getInitialGuess(void) const
Ptr< Vector< Real > > getEqualityMultiplier(void) const
Ptr< Constraint< Real > > getEqualityConstraint(void) const