75 Teuchos::RCP<const vector>
getVector(
const V& x ) {
76 using Teuchos::dyn_cast;
81 using Teuchos::dyn_cast;
92 return 100.0 * std::pow((*ex)[1] - std::pow((*ex)[0],2.0),2.0) + std::pow(1.0-(*ex)[0],2.0) +
93 90.0 * std::pow((*ex)[3] - std::pow((*ex)[2],2.0),2.0) + std::pow(1.0-(*ex)[2],2.0) +
94 10.1 * (std::pow((*ex)[1] - 1.0,2.0) + std::pow((*ex)[3]-1.0,2.0)) +
95 19.8 * ((*ex)[1] - 1.0) * ((*ex)[3] - 1.0);
104 (*eg)[0] = -4.0 * 100.0 * ((*ex)[1] - std::pow((*ex)[0],2.0)) * (*ex)[0] - 2.0 * (1.0-(*ex)[0]);
105 (*eg)[1] = 2.0 * 100.0 * ((*ex)[1] - std::pow((*ex)[0],2.0)) +
106 2.0 * 10.1 * ((*ex)[1] - 1.0) + 19.8*((*ex)[3] - 1.0);
107 (*eg)[2] = -4.0 * 90.0 * ((*ex)[3] - std::pow((*ex)[2],2.0)) * (*ex)[2] - 2.0 * (1.0-(*ex)[2]);
108 (*eg)[3] = 2.0 * 90.0 * ((*ex)[3] - std::pow((*ex)[2],2.0)) +
109 2.0 * 10.1 * ((*ex)[3] - 1.0) + 19.8*((*ex)[1] - 1.0);
119 Real h11 = -4.0 * 100.0 * (*ex)[1] + 12.0 * 100.0 * std::pow((*ex)[0],2.0) + 2.0;
120 Real h12 = -4.0 * 100.0 * (*ex)[0];
123 Real h21 = -4.0 * 100.0 * (*ex)[0];
124 Real h22 = 2.0 * 100.0 + 2.0 * 10.1;
129 Real h33 = -4.0 * 90.0 * (*ex)[3] + 12.0 * 90.0 * std::pow((*ex)[2],2.0) + 2.0;
130 Real h34 = -4.0 * 90.0 * (*ex)[2];
133 Real h43 = -4.0 * 90.0 * (*ex)[2];
134 Real h44 = 2.0 * 90.0 + 2.0 * 10.1;
136 (*ehv)[0] = h11 * (*ev)[0] + h12 * (*ev)[1] + h13 * (*ev)[2] + h14 * (*ev)[3];
137 (*ehv)[1] = h21 * (*ev)[0] + h22 * (*ev)[1] + h23 * (*ev)[2] + h24 * (*ev)[3];
138 (*ehv)[2] = h31 * (*ev)[0] + h32 * (*ev)[1] + h33 * (*ev)[2] + h34 * (*ev)[3];
139 (*ehv)[3] = h41 * (*ev)[0] + h42 * (*ev)[1] + h43 * (*ev)[2] + h44 * (*ev)[3];
148 typedef std::vector<Real>
vector;
153 using Teuchos::dyn_cast;
156 RCP<vector> x0p = dyn_cast<SV>(x0).
getVector();
157 RCP<vector> xp = dyn_cast<SV>(x).
getVector();
169 RCP<vector> l_rcp = rcp(
new vector(n,-10.0) );
170 RCP<vector> u_rcp = rcp(
new vector(n, 10.0) );
172 RCP<V> l = rcp(
new SV(l_rcp) );
173 RCP<V> u = rcp(
new SV(u_rcp) );
Provides the interface to evaluate objective functions.
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
Teuchos::RCP< const vector > getVector(const V &x)
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Contains definitions of custom data types in ROL.
std::vector< Real > vector
Defines the linear algebra or vector space interface.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
W. Hock and K. Schittkowski 38th test function.
void getHS38(Teuchos::RCP< Objective< Real > > &obj, Teuchos::RCP< BoundConstraint< Real > > &con, Vector< Real > &x0, Vector< Real > &x)
Provides the std::vector implementation of the ROL::Vector interface.
Teuchos::RCP< vector > getVector(V &x)
Provides the interface to apply upper and lower bound constraints.