49 #ifndef ROL_ZAKHAROV_HPP
50 #define ROL_ZAKHAROV_HPP
64 ROL::Ptr<Vector<Real> >
k_;
73 Real xdotx = x.
dot(x);
74 Real kdotx = x.
dot(*
k_);
76 Real val = xdotx + pow(kdotx,2)/4.0 + pow(kdotx,4)/16.0;
83 Real kdotx = x.
dot(*
k_);
84 Real coeff = 0.25*(2.0*kdotx+pow(kdotx,3.0));
93 Real kdotd = d.
dot(*
k_);
94 Real kdotx = x.
dot(*
k_);
95 Real xdotd = x.
dot(d);
97 Real coeff = 0.25*(2.0*kdotx+pow(kdotx,3.0));
99 Real deriv = 2*xdotd + coeff*kdotd;
108 Real kdotx = x.
dot(*
k_);
109 Real kdotv = v.
dot(*
k_);
110 Real coeff = 0.25*(2.0+3.0*pow(kdotx,2.0))*kdotv;
119 Real kdotv = v.
dot(*
k_);
120 Real kdotx = x.
dot(*
k_);
121 Real kdotk = (*k_).dot(*
k_);
122 Real coeff = -kdotv/(2.0*kdotk+16.0/(2.0+3.0*pow(kdotx,2.0)));
141 ROL::Ptr<std::vector<Real> > k_ptr = ROL::makePtr<std::vector<Real>>(n,0.0);
142 for (
int i = 0; i < n; i++ ) {
145 ROL::Ptr<Vector<Real> > k = ROL::makePtr<StdVector<Real>>(k_ptr);
146 return ROL::makePtr<Objective_Zakharov<Real>>(k);
153 ROL::Ptr<std::vector<Real> > x0p = ROL::makePtr<std::vector<Real>>(n,3.0);
154 return ROL::makePtr<StdVector<Real>>(x0p);
161 ROL::Ptr<std::vector<Real> > xp = ROL::makePtr<std::vector<Real>>(n,0.0);
162 return ROL::makePtr<StdVector<Real>>(xp);
Provides the interface to evaluate objective functions.
Objective_Zakharov(const ROL::Ptr< Vector< Real > > k)
virtual void scale(const Real alpha)=0
Compute where .
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
void invHessVec(Vector< Real > &ihv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
Ptr< Vector< Real > > getSolution(const int i=0) const
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Defines the linear algebra or vector space interface.
ROL::Ptr< Vector< Real > > k_
virtual Real dot(const Vector &x) const =0
Compute where .
Real value(const Vector< Real > &x, Real &tol)
Compute value.
Ptr< Objective< Real > > getObjective(void) const
Ptr< Vector< Real > > getInitialGuess(void) const
Contains definitions of test objective functions.
Real dirDeriv(const Vector< Real > &x, const Vector< Real > &d, Real &tol)
Compute directional derivative.
virtual void set(const Vector &x)
Set where .