ROL
ROL_Objective.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_OBJECTIVE_H
11#define ROL_OBJECTIVE_H
12
13#include "ROL_Vector.hpp"
14#include "ROL_UpdateType.hpp"
15#include "ROL_Types.hpp"
16#include <iostream>
17
41namespace ROL {
42
43template<typename Real>
44class Objective {
45private:
46 // Vector storage used for FD approximations (default are null pointers)
47 Ptr<Vector<Real>> prim_, dual_, basis_;
48
49public:
50
51 virtual ~Objective() {}
52
53 Objective() : prim_(nullPtr), dual_(nullPtr), basis_(nullPtr) {}
54
62 virtual void update( const Vector<Real> &x, UpdateType type, int iter = -1 ) {
63 ROL_UNUSED(x);
64 ROL_UNUSED(type);
65 ROL_UNUSED(iter);
66 }
67
75 virtual void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
76 ROL_UNUSED(x);
77 ROL_UNUSED(flag);
78 ROL_UNUSED(iter);
79 }
80
87 virtual Real value( const Vector<Real> &x, Real &tol ) = 0;
88
102 virtual void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) ;
103
111 virtual Real dirDeriv( const Vector<Real> &x, const Vector<Real> &d, Real &tol ) ;
112
121 virtual void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol );
122
131 virtual void invHessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
132 ROL_UNUSED(hv);
133 ROL_UNUSED(v);
134 ROL_UNUSED(x);
135 ROL_UNUSED(tol);
136 ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
137 ">>> ERROR (ROL::Objective): invHessVec not implemented!");
138 //hv.set(v.dual());
139 }
140
149 virtual void precond( Vector<Real> &Pv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
150 ROL_UNUSED(x);
151 ROL_UNUSED(tol);
152 Pv.set(v.dual());
153 }
154
163 virtual void prox( Vector<Real> &Pv, const Vector<Real> &v, Real t, Real &tol){
164 ROL_UNUSED(Pv);
165 ROL_UNUSED(v);
166 ROL_UNUSED(t);
167 ROL_UNUSED(tol);
168 ROL_TEST_FOR_EXCEPTION(true, std::invalid_argument,
169 ">>> ERROR (ROL::Objective): prox not implemented!");
170 }
171
172
182 virtual void proxJacVec( Vector<Real> &Jv, const Vector<Real> &v, const Vector<Real> &x, Real t, Real &tol);
183
204 virtual std::vector<std::vector<Real>> checkGradient( const Vector<Real> &x,
205 const Vector<Real> &d,
206 const bool printToStream = true,
207 std::ostream & outStream = std::cout,
208 const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
209 const int order = 1 ) {
210 return checkGradient(x, x.dual(), d, printToStream, outStream, numSteps, order);
211 }
212
235 virtual std::vector<std::vector<Real>> checkGradient( const Vector<Real> &x,
236 const Vector<Real> &g,
237 const Vector<Real> &d,
238 const bool printToStream = true,
239 std::ostream & outStream = std::cout,
240 const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
241 const int order = 1 );
242
243
264 virtual std::vector<std::vector<Real>> checkGradient( const Vector<Real> &x,
265 const Vector<Real> &d,
266 const std::vector<Real> &steps,
267 const bool printToStream = true,
268 std::ostream & outStream = std::cout,
269 const int order = 1 ) {
270
271 return checkGradient(x, x.dual(), d, steps, printToStream, outStream, order);
272
273 }
274
275
298 virtual std::vector<std::vector<Real>> checkGradient( const Vector<Real> &x,
299 const Vector<Real> &g,
300 const Vector<Real> &d,
301 const std::vector<Real> &steps,
302 const bool printToStream = true,
303 std::ostream & outStream = std::cout,
304 const int order = 1 );
305
326 virtual std::vector<std::vector<Real>> checkHessVec( const Vector<Real> &x,
327 const Vector<Real> &v,
328 const bool printToStream = true,
329 std::ostream & outStream = std::cout,
330 const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
331 const int order = 1 ) {
332
333 return checkHessVec(x, x.dual(), v, printToStream, outStream, numSteps, order);
334
335 }
336
358 virtual std::vector<std::vector<Real>> checkHessVec( const Vector<Real> &x,
359 const Vector<Real> &hv,
360 const Vector<Real> &v,
361 const bool printToStream = true,
362 std::ostream & outStream = std::cout,
363 const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
364 const int order = 1) ;
365
366
387 virtual std::vector<std::vector<Real>> checkHessVec( const Vector<Real> &x,
388 const Vector<Real> &v,
389 const std::vector<Real> &steps,
390 const bool printToStream = true,
391 std::ostream & outStream = std::cout,
392 const int order = 1 ) {
393
394 return checkHessVec(x, x.dual(), v, steps, printToStream, outStream, order);
395
396 }
397
419 virtual std::vector<std::vector<Real>> checkHessVec( const Vector<Real> &x,
420 const Vector<Real> &hv,
421 const Vector<Real> &v,
422 const std::vector<Real> &steps,
423 const bool printToStream = true,
424 std::ostream & outStream = std::cout,
425 const int order = 1) ;
426
427
442 virtual std::vector<Real> checkHessSym( const Vector<Real> &x,
443 const Vector<Real> &v,
444 const Vector<Real> &w,
445 const bool printToStream = true,
446 std::ostream & outStream = std::cout ) {
447
448 return checkHessSym(x, x.dual(), v, w, printToStream, outStream);
449
450 }
451
467 virtual std::vector<Real> checkHessSym( const Vector<Real> &x,
468 const Vector<Real> &hv,
469 const Vector<Real> &v,
470 const Vector<Real> &w,
471 const bool printToStream = true,
472 std::ostream & outStream = std::cout );
473
492 virtual std::vector<std::vector<Real>> checkProxJacVec( const Vector<Real> &x,
493 const Vector<Real> &v,
494 Real t = Real(1),
495 bool printToStream = true,
496 std::ostream &outStream = std::cout,
497 int numSteps = ROL_NUM_CHECKDERIV_STEPS);
498
499// Definitions for parametrized (stochastic) objective functions
500private:
501 std::vector<Real> param_;
502
503protected:
504 const std::vector<Real> getParameter(void) const {
505 return param_;
506 }
507
508public:
509 virtual void setParameter(const std::vector<Real> &param) {
510 param_.assign(param.begin(),param.end());
511 }
512
513}; // class Objective
514
515} // namespace ROL
516
517// include templated definitions
518#include <ROL_ObjectiveDef.hpp>
519
520#endif
Contains definitions of custom data types in ROL.
#define ROL_NUM_CHECKDERIV_STEPS
Number of steps for derivative checks.
Definition ROL_Types.hpp:40
#define ROL_UNUSED(x)
Provides the interface to evaluate objective functions.
virtual std::vector< Real > checkHessSym(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout)
Hessian symmetry check.
virtual Real dirDeriv(const Vector< Real > &x, const Vector< Real > &d, Real &tol)
Compute directional derivative.
virtual void setParameter(const std::vector< Real > &param)
Ptr< Vector< Real > > prim_
virtual void prox(Vector< Real > &Pv, const Vector< Real > &v, Real t, Real &tol)
Compute the proximity operator.
virtual std::vector< std::vector< Real > > checkGradient(const Vector< Real > &x, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference gradient check.
const std::vector< Real > getParameter(void) const
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
std::vector< Real > param_
virtual ~Objective()
virtual std::vector< std::vector< Real > > checkProxJacVec(const Vector< Real > &x, const Vector< Real > &v, Real t=Real(1), bool printToStream=true, std::ostream &outStream=std::cout, int numSteps=ROL_NUM_CHECKDERIV_STEPS)
Finite-difference proximity operator Jacobian-applied-to-vector check.
virtual std::vector< std::vector< Real > > checkHessVec(const Vector< Real > &x, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference Hessian-applied-to-vector check with specified step sizes.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
Ptr< Vector< Real > > dual_
virtual void proxJacVec(Vector< Real > &Jv, const Vector< Real > &v, const Vector< Real > &x, Real t, Real &tol)
Apply the Jacobian of the proximity operator.
virtual std::vector< std::vector< Real > > checkGradient(const Vector< Real > &x, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference gradient check with specified step sizes.
virtual std::vector< std::vector< Real > > checkHessVec(const Vector< Real > &x, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference Hessian-applied-to-vector check.
Ptr< Vector< Real > > basis_
virtual void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply preconditioner to vector.
Defines the linear algebra or vector space interface.
virtual void set(const Vector &x)
Set where .
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...