ROL
ROL_ParaboloidCircle.hpp
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
18
19#ifndef ROL_PARABOLOIDCIRCLE_HPP
20#define ROL_PARABOLOIDCIRCLE_HPP
21
22#include "ROL_TestProblem.hpp"
23#include "ROL_StdVector.hpp"
24#include "Teuchos_SerialDenseVector.hpp"
25#include "Teuchos_SerialDenseSolver.hpp"
26
27namespace ROL {
28namespace ZOO {
29
33 template< class Real, class XPrim=StdVector<Real>, class XDual=StdVector<Real> >
35
36 typedef std::vector<Real> vector;
37 typedef Vector<Real> V;
38
39 typedef typename vector::size_type uint;
40
41
42 private:
43
44 template<class VectorType>
45 ROL::Ptr<const vector> getVector( const V& x ) {
46
47 return dynamic_cast<const VectorType&>(x).getVector();
48 }
49
50 template<class VectorType>
51 ROL::Ptr<vector> getVector( V& x ) {
52
53 return dynamic_cast<VectorType&>(x).getVector();
54 }
55
56 public:
58
59 Real value( const Vector<Real> &x, Real &tol ) {
60
61
62 ROL::Ptr<const vector> xp = getVector<XPrim>(x);
63
64 uint n = xp->size();
65 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, objective value): "
66 "Primal vector x must be of length 2.");
67
68 Real x1 = (*xp)[0];
69 Real x2 = (*xp)[1];
70
71 Real val = x1*x1 + x2*x2;
72
73 return val;
74 }
75
76 void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
77
78
79 ROL::Ptr<const vector> xp = getVector<XPrim>(x);
80 ROL::Ptr<vector> gp = getVector<XDual>(g);
81
82 uint n = xp->size();
83 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, objective gradient): "
84 " Primal vector x must be of length 2.");
85
86 n = gp->size();
87 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, objective gradient): "
88 "Gradient vector g must be of length 2.");
89
90 Real x1 = (*xp)[0];
91 Real x2 = (*xp)[1];
92
93 Real two(2);
94
95 (*gp)[0] = two*x1;
96 (*gp)[1] = two*x2;
97 }
98
99 void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
100
101
102 ROL::Ptr<const vector> xp = getVector<XPrim>(x);
103 ROL::Ptr<const vector> vp = getVector<XPrim>(v);
104 ROL::Ptr<vector> hvp = getVector<XDual>(hv);
105
106 uint n = xp->size();
107 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, objective hessVec): "
108 "Primal vector x must be of length 2.");
109
110 n = vp->size();
111 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, objective hessVec): "
112 "Input vector v must be of length 2.");
113
114 n = hvp->size();
115 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, objective hessVec): "
116 "Output vector hv must be of length 2.");
117
118 Real v1 = (*vp)[0];
119 Real v2 = (*vp)[1];
120
121 Real two(2);
122
123 (*hvp)[0] = two*v1;
124 (*hvp)[1] = two*v2;
125 }
126
127 };
128
129
132 template<class Real, class XPrim=StdVector<Real>, class XDual=StdVector<Real>, class CPrim=StdVector<Real>, class CDual=StdVector<Real> >
134
135 typedef std::vector<Real> vector;
137
138 typedef typename vector::size_type uint;
139
140 private:
141 template<class VectorType>
142 ROL::Ptr<const vector> getVector( const V& x ) {
143
144 return dynamic_cast<const VectorType&>(x).getVector();
145 }
146
147 template<class VectorType>
148 ROL::Ptr<vector> getVector( V& x ) {
149
150 return dynamic_cast<VectorType&>(x).getVector();
151 }
152
153 public:
155
156 void value( Vector<Real> &c, const Vector<Real> &x, Real &tol ) {
157
158
159 ROL::Ptr<const vector> xp = getVector<XPrim>(x);
160 ROL::Ptr<vector> cp = getVector<CPrim>(c);
161
162 uint n = xp->size();
163 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint value): "
164 "Primal vector x must be of length 2.");
165
166 uint m = cp->size();
167 ROL_TEST_FOR_EXCEPTION( (m != 1), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint value): "
168 "Constraint vector c must be of length 1.");
169
170 Real x1 = (*xp)[0];
171 Real x2 = (*xp)[1];
172
173 Real one(1), two(2);
174
175 (*cp)[0] = (x1-two)*(x1-two) + x2*x2 - one;
176 }
177
178 void applyJacobian( Vector<Real> &jv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
179
180
181 ROL::Ptr<const vector> xp = getVector<XPrim>(x);
182 ROL::Ptr<const vector> vp = getVector<XPrim>(v);
183 ROL::Ptr<vector> jvp = getVector<CPrim>(jv);
184
185 uint n = xp->size();
186 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint applyJacobian): "
187 "Primal vector x must be of length 2.");
188
189 uint d = vp->size();
190 ROL_TEST_FOR_EXCEPTION( (d != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint applyJacobian): "
191 "Input vector v must be of length 2.");
192 d = jvp->size();
193 ROL_TEST_FOR_EXCEPTION( (d != 1), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint applyJacobian): "
194 "Output vector jv must be of length 1.");
195
196 Real x1 = (*xp)[0];
197 Real x2 = (*xp)[1];
198
199 Real v1 = (*vp)[0];
200 Real v2 = (*vp)[1];
201
202 Real two(2);
203
204 (*jvp)[0] = two*(x1-two)*v1 + two*x2*v2;
205 } //applyJacobian
206
207 void applyAdjointJacobian( Vector<Real> &ajv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
208
209
210 ROL::Ptr<const vector> xp = getVector<XPrim>(x);
211 ROL::Ptr<const vector> vp = getVector<CDual>(v);
212 ROL::Ptr<vector> ajvp = getVector<XDual>(ajv);
213
214 uint n = xp->size();
215 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointJacobian): "
216 "Primal vector x must be of length 2.");
217
218 uint d = vp->size();
219 ROL_TEST_FOR_EXCEPTION( (d != 1), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointJacobian): "
220 "Input vector v must be of length 1.");
221
222 d = ajvp->size();
223 ROL_TEST_FOR_EXCEPTION( (d != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointJacobian): "
224 "Output vector ajv must be of length 2.");
225
226 Real x1 = (*xp)[0];
227 Real x2 = (*xp)[1];
228
229 Real v1 = (*vp)[0];
230
231 Real two(2);
232
233 (*ajvp)[0] = two*(x1-two)*v1;
234 (*ajvp)[1] = two*x2*v1;
235
236 } //applyAdjointJacobian
237
238 void applyAdjointHessian( Vector<Real> &ahuv, const Vector<Real> &u, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
239
240 bool useFD = true;
241
242 if (useFD) {
243 Constraint<Real>::applyAdjointHessian( ahuv, u, v, x, tol );
244 }
245 else {
246
247 ROL::Ptr<const vector> xp = getVector<XPrim>(x);
248 ROL::Ptr<const vector> up = getVector<CDual>(u);
249 ROL::Ptr<const vector> vp = getVector<XPrim>(v);
250 ROL::Ptr<vector> ahuvp = getVector<XDual>(ahuv);
251
252 uint n = xp->size();
253 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
254 "Primal vector x must be of length 2.");
255
256 n = vp->size();
257 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
258 "Direction vector v must be of length 2.");
259
260 n = ahuvp->size();
261 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
262 "Output vector ahuv must be of length 2.");
263 uint d = up->size();
264 ROL_TEST_FOR_EXCEPTION( (d != 1), std::invalid_argument, ">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
265 "Dual constraint vector u must be of length 1.");
266
267 Real v1 = (*vp)[0];
268 Real v2 = (*vp)[1];
269
270 Real u1 = (*up)[0];
271
272 Real two(2);
273
274 (*ahuvp)[0] = two*u1*v1;
275 (*ahuvp)[1] = two*u1*v2;
276 }
277 } //applyAdjointHessian
278
279 };
280
281
282 template<class Real, class XPrim=StdVector<Real>, class XDual=StdVector<Real>, class CPrim=StdVector<Real>, class CDual=StdVector<Real> >
283 class getParaboloidCircle : public TestProblem<Real> {
284 typedef std::vector<Real> vector;
285 typedef typename vector::size_type uint;
286 public:
288
289 Ptr<Objective<Real>> getObjective(void) const {
290 // Instantiate objective function.
291 return ROL::makePtr<Objective_ParaboloidCircle<Real,XPrim,XDual>>();
292 }
293
294 Ptr<Vector<Real>> getInitialGuess(void) const {
295 uint n = 2;
296 // Get initial guess.
297 ROL::Ptr<vector> x0p = makePtr<vector>(n,0.0);
298 (*x0p)[0] = static_cast<Real>(rand())/static_cast<Real>(RAND_MAX);
299 (*x0p)[1] = static_cast<Real>(rand())/static_cast<Real>(RAND_MAX);
300 return makePtr<XPrim>(x0p);
301 }
302
303 Ptr<Vector<Real>> getSolution(const int i = 0) const {
304 uint n = 2;
305 // Get solution.
306 Real zero(0), one(1);
307 ROL::Ptr<vector> solp = makePtr<vector>(n,0.0);
308 (*solp)[0] = one;
309 (*solp)[1] = zero;
310 return makePtr<XPrim>(solp);
311 }
312
313 Ptr<Constraint<Real>> getEqualityConstraint(void) const {
314 // Instantiate constraints.
315 return ROL::makePtr<Constraint_ParaboloidCircle<Real,XPrim,XDual,CPrim,CDual>>();
316 }
317
318 Ptr<Vector<Real>> getEqualityMultiplier(void) const {
319 ROL::Ptr<vector> lp = makePtr<vector>(1,0.0);
320 return makePtr<CDual>(lp);
321 }
322 };
323
324} // End ZOO Namespace
325} // End ROL Namespace
326
327#endif
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0 zero)()
Contains definitions of test objective functions.
virtual void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the derivative of the adjoint of the constraint Jacobian at to vector in direction ,...
Defines the linear algebra or vector space interface.
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
ROL::Ptr< const vector > getVector(const V &x)
void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Real value(const Vector< Real > &x, Real &tol)
ROL::Ptr< const vector > getVector(const V &x)
Ptr< Vector< Real > > getSolution(const int i=0) const
Ptr< Constraint< Real > > getEqualityConstraint(void) const
Ptr< Objective< Real > > getObjective(void) const
Ptr< Vector< Real > > getEqualityMultiplier(void) const
Ptr< Vector< Real > > getInitialGuess(void) const