53#ifndef ROL_PARABOLOIDCIRCLE_HPP
54#define ROL_PARABOLOIDCIRCLE_HPP
58#include "Teuchos_SerialDenseVector.hpp"
59#include "Teuchos_SerialDenseSolver.hpp"
67 template<
class Real,
class XPrim=StdVector<Real>,
class XDual=StdVector<Real> >
73 typedef typename vector::size_type
uint;
78 template<
class VectorType>
81 return dynamic_cast<const VectorType&
>(x).
getVector();
84 template<
class VectorType>
87 return dynamic_cast<VectorType&
>(x).
getVector();
99 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective value): "
100 "Primal vector x must be of length 2.");
105 Real val = x1*x1 + x2*x2;
117 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective gradient): "
118 " Primal vector x must be of length 2.");
121 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective gradient): "
122 "Gradient vector g must be of length 2.");
141 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective hessVec): "
142 "Primal vector x must be of length 2.");
145 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective hessVec): "
146 "Input vector v must be of length 2.");
149 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective hessVec): "
150 "Output vector hv must be of length 2.");
166 template<
class Real,
class XPrim=StdVector<Real>,
class XDual=StdVector<Real>,
class CPrim=StdVector<Real>,
class CDual=StdVector<Real> >
172 typedef typename vector::size_type
uint;
175 template<
class VectorType>
178 return dynamic_cast<const VectorType&
>(x).
getVector();
181 template<
class VectorType>
184 return dynamic_cast<VectorType&
>(x).
getVector();
197 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint value): "
198 "Primal vector x must be of length 2.");
201 ROL_TEST_FOR_EXCEPTION( (m != 1), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint value): "
202 "Constraint vector c must be of length 1.");
209 (*cp)[0] = (x1-two)*(x1-two) + x2*x2 - one;
220 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyJacobian): "
221 "Primal vector x must be of length 2.");
224 ROL_TEST_FOR_EXCEPTION( (d != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyJacobian): "
225 "Input vector v must be of length 2.");
227 ROL_TEST_FOR_EXCEPTION( (d != 1), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyJacobian): "
228 "Output vector jv must be of length 1.");
238 (*jvp)[0] = two*(x1-two)*v1 + two*x2*v2;
249 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointJacobian): "
250 "Primal vector x must be of length 2.");
253 ROL_TEST_FOR_EXCEPTION( (d != 1), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointJacobian): "
254 "Input vector v must be of length 1.");
257 ROL_TEST_FOR_EXCEPTION( (d != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointJacobian): "
258 "Output vector ajv must be of length 2.");
267 (*ajvp)[0] = two*(x1-two)*v1;
268 (*ajvp)[1] = two*x2*v1;
287 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
288 "Primal vector x must be of length 2.");
291 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
292 "Direction vector v must be of length 2.");
295 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
296 "Output vector ahuv must be of length 2.");
298 ROL_TEST_FOR_EXCEPTION( (d != 1), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
299 "Dual constraint vector u must be of length 1.");
308 (*ahuvp)[0] = two*u1*v1;
309 (*ahuvp)[1] = two*u1*v2;
316 template<
class Real,
class XPrim=StdVector<Real>,
class XDual=StdVector<Real>,
class CPrim=StdVector<Real>,
class CDual=StdVector<Real> >
319 typedef typename vector::size_type
uint;
325 return ROL::makePtr<Objective_ParaboloidCircle<Real,XPrim,XDual>>();
331 ROL::Ptr<vector> x0p = makePtr<vector>(n,0.0);
332 (*x0p)[0] =
static_cast<Real
>(rand())/
static_cast<Real
>(RAND_MAX);
333 (*x0p)[1] =
static_cast<Real
>(rand())/
static_cast<Real
>(RAND_MAX);
334 return makePtr<XPrim>(x0p);
340 Real
zero(0), one(1);
341 ROL::Ptr<vector> solp = makePtr<vector>(n,0.0);
344 return makePtr<XPrim>(solp);
349 return ROL::makePtr<Constraint_ParaboloidCircle<Real,XPrim,XDual,CPrim,CDual>>();
353 ROL::Ptr<vector> lp = makePtr<vector>(1,0.0);
354 return makePtr<CDual>(lp);
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0 zero)()
Contains definitions of test objective functions.
Defines the general constraint operator interface.
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 ,...
Provides the interface to evaluate objective functions.
Defines the linear algebra or vector space interface.
constraint c(x,y) = (x-2)^2 + y^2 - 1.
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Constraint_ParaboloidCircle()
ROL::Ptr< vector > getVector(V &x)
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
ROL::Ptr< const vector > getVector(const V &x)
std::vector< Real > vector
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)
Objective function: f(x,y) = x^2 + y^2.
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Objective_ParaboloidCircle()
std::vector< Real > vector
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
ROL::Ptr< vector > getVector(V &x)
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
getParaboloidCircle(void)
Ptr< Objective< Real > > getObjective(void) const
std::vector< Real > vector
Ptr< Vector< Real > > getEqualityMultiplier(void) const
Ptr< Vector< Real > > getInitialGuess(void) const