19#ifndef ROL_PARABOLOIDCIRCLE_HPP
20#define ROL_PARABOLOIDCIRCLE_HPP
24#include "Teuchos_SerialDenseVector.hpp"
25#include "Teuchos_SerialDenseSolver.hpp"
33 template<
class Real,
class XPrim=StdVector<Real>,
class XDual=StdVector<Real> >
39 typedef typename vector::size_type
uint;
44 template<
class VectorType>
47 return dynamic_cast<const VectorType&
>(x).
getVector();
50 template<
class VectorType>
53 return dynamic_cast<VectorType&
>(x).
getVector();
65 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective value): "
66 "Primal vector x must be of length 2.");
71 Real val = x1*x1 + x2*x2;
83 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective gradient): "
84 " Primal vector x must be of length 2.");
87 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective gradient): "
88 "Gradient vector g must be of length 2.");
107 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective hessVec): "
108 "Primal vector x must be of length 2.");
111 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective hessVec): "
112 "Input vector v must be of length 2.");
115 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, objective hessVec): "
116 "Output vector hv must be of length 2.");
132 template<
class Real,
class XPrim=StdVector<Real>,
class XDual=StdVector<Real>,
class CPrim=StdVector<Real>,
class CDual=StdVector<Real> >
138 typedef typename vector::size_type
uint;
141 template<
class VectorType>
144 return dynamic_cast<const VectorType&
>(x).
getVector();
147 template<
class VectorType>
150 return dynamic_cast<VectorType&
>(x).
getVector();
163 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint value): "
164 "Primal vector x must be of length 2.");
167 ROL_TEST_FOR_EXCEPTION( (m != 1), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint value): "
168 "Constraint vector c must be of length 1.");
175 (*cp)[0] = (x1-two)*(x1-two) + x2*x2 - one;
186 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyJacobian): "
187 "Primal vector x must be of length 2.");
190 ROL_TEST_FOR_EXCEPTION( (d != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyJacobian): "
191 "Input vector v must be of length 2.");
193 ROL_TEST_FOR_EXCEPTION( (d != 1), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyJacobian): "
194 "Output vector jv must be of length 1.");
204 (*jvp)[0] = two*(x1-two)*v1 + two*x2*v2;
215 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointJacobian): "
216 "Primal vector x must be of length 2.");
219 ROL_TEST_FOR_EXCEPTION( (d != 1), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointJacobian): "
220 "Input vector v must be of length 1.");
223 ROL_TEST_FOR_EXCEPTION( (d != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointJacobian): "
224 "Output vector ajv must be of length 2.");
233 (*ajvp)[0] = two*(x1-two)*v1;
234 (*ajvp)[1] = two*x2*v1;
253 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
254 "Primal vector x must be of length 2.");
257 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
258 "Direction vector v must be of length 2.");
261 ROL_TEST_FOR_EXCEPTION( (n != 2), std::invalid_argument,
">>> ERROR (ROL_ParaboloidCircle, constraint applyAdjointHessian): "
262 "Output vector ahuv must be of length 2.");
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.");
274 (*ahuvp)[0] = two*u1*v1;
275 (*ahuvp)[1] = two*u1*v2;
282 template<
class Real,
class XPrim=StdVector<Real>,
class XDual=StdVector<Real>,
class CPrim=StdVector<Real>,
class CDual=StdVector<Real> >
285 typedef typename vector::size_type
uint;
291 return ROL::makePtr<Objective_ParaboloidCircle<Real,XPrim,XDual>>();
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);
306 Real
zero(0), one(1);
307 ROL::Ptr<vector> solp = makePtr<vector>(n,0.0);
310 return makePtr<XPrim>(solp);
315 return ROL::makePtr<Constraint_ParaboloidCircle<Real,XPrim,XDual,CPrim,CDual>>();
319 ROL::Ptr<vector> lp = makePtr<vector>(1,0.0);
320 return makePtr<CDual>(lp);
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)
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)
void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
std::vector< Real > vector
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)
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
Ptr< Vector< Real > > getEqualityMultiplier(void) const
Ptr< Vector< Real > > getInitialGuess(void) const
std::vector< Real > vector