80     ROL::Ptr<const vector> xp = 
dynamic_cast<const SV&
>(x).getVector();
 
   85     ROL::Ptr<const vector> xp = 
dynamic_cast<const SV&
>(x).getVector();
 
   86     ROL::Ptr<vector> gp = 
dynamic_cast<SV&
>(g).getVector();
 
  112     ROL::Ptr<vector> cp = 
dynamic_cast<SV&
>(c).getVector();
 
  113     ROL::Ptr<const vector> xp = 
dynamic_cast<const SV&
>(x).getVector();
 
  115     (*cp)[0] = (*xp)[1]-std::pow((*xp)[0],3)-std::pow((*xp)[2],2);
 
  121     ROL::Ptr<vector> jvp = 
dynamic_cast<SV&
>(jv).getVector();
 
  122     ROL::Ptr<const vector> vp = 
dynamic_cast<const SV&
>(v).getVector();
 
  123     ROL::Ptr<const vector> xp = 
dynamic_cast<const SV&
>(x).getVector();
 
  125     (*jvp)[0] = (*vp)[1] - 3.0*(*xp)[0]*(*xp)[0]*(*vp)[0] - 2.0*(*xp)[2]*(*vp)[2];
 
  132     ROL::Ptr<vector> ajvp = 
dynamic_cast<SV&
>(ajv).getVector();
 
  133     ROL::Ptr<const vector> vp = 
dynamic_cast<const SV&
>(v).getVector();
 
  134     ROL::Ptr<const vector> xp = 
dynamic_cast<const SV&
>(x).getVector();
 
  136     (*ajvp)[0] = -3.0*(*xp)[0]*(*xp)[0]*(*vp)[0];
 
  137     (*ajvp)[1] = (*vp)[0];
 
  138     (*ajvp)[2] = -2.0*(*xp)[2]*(*vp)[0];
 
  146     ROL::Ptr<vector> ahuvp = 
dynamic_cast<SV&
>(ahuv).getVector();
 
  147     ROL::Ptr<const vector> up = 
dynamic_cast<const SV&
>(u).getVector();
 
  148     ROL::Ptr<const vector> vp = 
dynamic_cast<const SV&
>(v).getVector();
 
  149     ROL::Ptr<const vector> xp = 
dynamic_cast<const SV&
>(x).getVector();
 
  151     (*ahuvp)[0] = -6.0*(*up)[0]*(*xp)[0]*(*vp)[0]; 
 
  153     (*ahuvp)[2] = -2.0*(*up)[0]*(*vp)[2];
 
  172     ROL::Ptr<vector> cp = 
dynamic_cast<SV&
>(c).getVector();
 
  173     ROL::Ptr<const vector> xp = 
dynamic_cast<const SV&
>(x).getVector();
 
  175     (*cp)[0] = std::pow((*xp)[0],2)-(*xp)[1]-std::pow((*xp)[3],2);
 
  181     ROL::Ptr<vector> jvp = 
dynamic_cast<SV&
>(jv).getVector();
 
  182     ROL::Ptr<const vector> vp = 
dynamic_cast<const SV&
>(v).getVector();
 
  183     ROL::Ptr<const vector> xp = 
dynamic_cast<const SV&
>(x).getVector();
 
  185     (*jvp)[0] = 2.0*(*xp)[0]*(*vp)[0] - (*vp)[1] - 2.0*(*xp)[3]*(*vp)[3];
 
  192     ROL::Ptr<vector> ajvp = 
dynamic_cast<SV&
>(ajv).getVector();
 
  193     ROL::Ptr<const vector> vp = 
dynamic_cast<const SV&
>(v).getVector();
 
  194     ROL::Ptr<const vector> xp = 
dynamic_cast<const SV&
>(x).getVector();
 
  196     (*ajvp)[0] = 2.0*(*xp)[0]*(*vp)[0];
 
  197     (*ajvp)[1] = -(*vp)[0]; 
 
  199     (*ajvp)[3] = -2.0*(*vp)[0]*(*xp)[3];
 
  206     ROL::Ptr<vector> ahuvp = 
dynamic_cast<SV&
>(ahuv).getVector();
 
  207     ROL::Ptr<const vector> up = 
dynamic_cast<const SV&
>(u).getVector();
 
  208     ROL::Ptr<const vector> vp = 
dynamic_cast<const SV&
>(v).getVector();
 
  209     ROL::Ptr<const vector> xp = 
dynamic_cast<const SV&
>(x).getVector();
 
  213     (*ahuvp)[0] = 2.0*(*up)[0]*(*vp)[0];
 
  216     (*ahuvp)[3] = -2.0*(*up)[0]*(*vp)[3];
 
  229     return ROL::makePtr<Objective_HS39<Real>>();
 
  236     ROL::Ptr<std::vector<Real> > x0p = ROL::makePtr<std::vector<Real>>(n,2.0);
 
  237     return ROL::makePtr<StdVector<Real>>(x0p);
 
  244     ROL::Ptr<std::vector<Real> > xp = ROL::makePtr<std::vector<Real>>(n,0.0);
 
  245     (*xp)[0] = 
static_cast<Real
>(1);
 
  246     (*xp)[1] = 
static_cast<Real
>(1);
 
  247     (*xp)[2] = 
static_cast<Real
>(0);
 
  248     (*xp)[3] = 
static_cast<Real
>(0);
 
  249     return ROL::makePtr<StdVector<Real>>(xp);
 
  253     std::vector<Ptr<Constraint<Real>>> cvec(2);
 
  254     cvec[0] = makePtr<Constraint_HS39a<Real>>();
 
  255     cvec[1] = makePtr<Constraint_HS39b<Real>>();
 
  256     return ROL::makePtr<Constraint_Partitioned<Real>>(cvec);
 
  260     std::vector<Ptr<Vector<Real>>> lvec(2);
 
  261     lvec[0] = makePtr<StdVector<Real>>(makePtr<std::vector<Real>>(1,0.0));
 
  262     lvec[1] = makePtr<StdVector<Real>>(makePtr<std::vector<Real>>(1,0.0));
 
  263     return ROL::makePtr<PartitionedVector<Real>>(lvec);
 
Provides the interface to evaluate objective functions. 
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector . 
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 ...
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector. 
void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the adjoint of the the constraint Jacobian at , , to vector . 
Ptr< Vector< Real > > getSolution(const int i=0) const 
Contains definitions of custom data types in ROL. 
Ptr< Objective< Real > > getObjective(void) const 
virtual void zero()
Set to zero vector. 
Defines the linear algebra or vector space interface. 
std::vector< Real > vector
std::vector< Real > vector
std::vector< Real > vector
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient. 
Ptr< Vector< Real > > getInitialGuess(void) const 
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 ...
Contains definitions of test objective functions. 
Ptr< Vector< Real > > getEqualityMultiplier(void) const 
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator  at . 
Ptr< Constraint< Real > > getEqualityConstraint(void) const 
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator  at . 
Real value(const Vector< Real > &x, Real &tol)
Compute value. 
void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the adjoint of the the constraint Jacobian at , , to vector . 
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector . 
Defines the general constraint operator interface. 
W. Hock and K. Schittkowski 39th test function.