44 #ifndef ROL_STDOBJECTIVE_H
45 #define ROL_STDOBJECTIVE_H
62 Real
sgn(
const Real x)
const {
63 const Real
zero(0), one(1);
64 return (x <
zero ? -one : one);
68 virtual void update(
const std::vector<Real> &x,
bool flag =
true,
int iter = -1 ) {}
76 virtual Real
value(
const std::vector<Real> &x, Real &tol ) = 0;
84 virtual void gradient( std::vector<Real> &g,
const std::vector<Real> &x, Real &tol ) {
85 const unsigned size = x.size();
86 std::vector<Real> y; y.assign(x.begin(),x.end());
87 const Real cbrteps = std::cbrt(ROL::ROL_EPSILON<Real>()), one(1);
89 const Real val =
value(x,tol);
90 for (
unsigned i = 0; i < size; ++i) {
92 h = cbrteps * std::max(std::abs(xi),one) *
sgn(xi);
96 g[i] = (
value(y,tol) - val)/h;
109 virtual Real
dirDeriv(
const std::vector<Real> &x,
const std::vector<Real> &d, Real &tol ) {
110 ROL_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
111 ">>> ERROR (ROL::StdObjective): dirDeriv not implemented!");
121 catch (std::exception &e) {
126 virtual void hessVec( std::vector<Real> &hv,
const std::vector<Real> &v,
const std::vector<Real> &x, Real &tol ) {
127 ROL_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
128 ">>> ERROR (ROL::StdObjective): hessVec not implemented!");
139 catch (std::exception &e) {
144 virtual void invHessVec( std::vector<Real> &hv,
const std::vector<Real> &v,
const std::vector<Real> &x, Real &tol ) {
145 ROL_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
146 ">>> ERROR (ROL::StdObjective): invHessVec not implemented!");
157 virtual void precond( std::vector<Real> &Pv,
const std::vector<Real> &v,
const std::vector<Real> &x, Real &tol ) {
158 Pv.assign(v.begin(),v.end());
virtual Real value(const std::vector< Real > &x, Real &tol)=0
Provides the interface to evaluate objective functions.
Real sgn(const Real x) const
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
virtual Real dirDeriv(const Vector< Real > &x, const Vector< Real > &d, Real &tol)
Compute directional derivative.
virtual void hessVec(std::vector< Real > &hv, const std::vector< Real > &v, const std::vector< Real > &x, Real &tol)
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Ptr< const std::vector< Element > > getVector() const
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
virtual Real dirDeriv(const std::vector< Real > &x, const std::vector< Real > &d, Real &tol)
Defines the linear algebra or vector space interface.
Real dirDeriv(const Vector< Real > &x, const Vector< Real > &d, Real &tol)
Compute directional derivative.
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Real value(const Vector< Real > &x, Real &tol)
Compute value.
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
Specializes the ROL::Objective interface for objective functions that operate on ROL::StdVector's.
void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply preconditioner to vector.
virtual void precond(std::vector< Real > &Pv, const std::vector< Real > &v, const std::vector< Real > &x, Real &tol)
virtual void invHessVec(std::vector< Real > &hv, const std::vector< Real > &v, const std::vector< Real > &x, Real &tol)
virtual void update(const std::vector< Real > &x, bool flag=true, int iter=-1)
virtual void gradient(std::vector< Real > &g, const std::vector< Real > &x, Real &tol)