ROL
ROL::Objective< Real > Class Template Referenceabstract

Provides the interface to evaluate objective functions. More...

#include <ROL_Objective.hpp>

Inheritance diagram for ROL::Objective< Real >:

Public Member Functions

virtual ~Objective ()
 Objective ()
virtual void update (const Vector< Real > &x, UpdateType type, int iter=-1)
 Update objective function.
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 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 (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply Hessian approximation to vector.
virtual void invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply inverse Hessian approximation to vector.
virtual void precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply preconditioner to vector.
virtual void prox (Vector< Real > &Pv, const Vector< Real > &v, Real t, Real &tol)
 Compute the proximity operator.
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 bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference gradient check.
virtual std::vector< std::vector< Real > > checkGradient (const Vector< Real > &x, const Vector< Real > &g, 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.
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 > > checkGradient (const Vector< Real > &x, const Vector< Real > &g, 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.
virtual std::vector< std::vector< Real > > checkHessVec (const Vector< Real > &x, const Vector< Real > &hv, 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.
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 std::vector< std::vector< Real > > checkHessVec (const Vector< Real > &x, const Vector< Real > &hv, 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 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 std::vector< Real > checkHessSym (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout)
 Hessian symmetry check.
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 void setParameter (const std::vector< Real > &param)

Protected Member Functions

const std::vector< Real > getParameter (void) const

Private Attributes

Ptr< Vector< Real > > prim_
Ptr< Vector< Real > > dual_
Ptr< Vector< Real > > basis_
std::vector< Real > param_

Detailed Description

template<typename Real>
class ROL::Objective< Real >

Provides the interface to evaluate objective functions.

Provides the definition of the objective function interface.

ROL's objective function interface is designed for Fr$eacute;chet differentiable functionals \(f:\mathcal{X}\to\mathbb{R}\), where \(\mathcal{X}\) is a Banach space. The basic operator interace, to be implemented by the user, requires:

  • value – objective function evaluation.

It is strongly recommended that the user additionally overload:

  • gradient – the objective function gradient – the default is a finite-difference approximation;
  • hessVec – the action of the Hessian – the default is a finite-difference approximation.

The user may also overload:

  • update – update the objective function at each new iteration;
  • dirDeriv – compute the directional derivative – the default is a finite-difference approximation;
  • invHessVec – the action of the inverse Hessian;
  • precond – the action of a preconditioner for the Hessian.

Definition at line 44 of file ROL_Objective.hpp.

Constructor & Destructor Documentation

◆ ~Objective()

template<typename Real>
virtual ROL::Objective< Real >::~Objective ( )
inlinevirtual

Definition at line 51 of file ROL_Objective.hpp.

◆ Objective()

template<typename Real>
ROL::Objective< Real >::Objective ( )
inline

Definition at line 53 of file ROL_Objective.hpp.

References basis_, dual_, and prim_.

Member Function Documentation

◆ update() [1/2]

template<typename Real>
virtual void ROL::Objective< Real >::update ( const Vector< Real > & x,
UpdateType type,
int iter = -1 )
inlinevirtual

Update objective function.

This function updates the objective function at new iterations.

Parameters
[in]xis the new iterate.
[in]typeis the type of update requested.
[in]iteris the outer algorithm iterations count.

Definition at line 62 of file ROL_Objective.hpp.

References ROL_UNUSED.

Referenced by ROL::CompositeStep< Real >::accept(), ROL::BundleStep< Real >::compute(), ROL::TypeB::LinMoreAlgorithm< Real >::computeValue(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::computeValue(), ROL::TypeP::TrustRegionAlgorithm< Real >::computeValue(), ROL::TypeP::TrustRegionAlgorithm< Real >::dbls(), ROL::TypeP::TrustRegionAlgorithm< Real >::dcauchy(), ROL::LineSearch_U< Real >::dirDeriv(), ROL::TypeP::TrustRegionAlgorithm< Real >::dncg(), ROL::TypeP::TrustRegionAlgorithm< Real >::dspg(), ROL::TypeP::TrustRegionAlgorithm< Real >::dspg2(), ROL::TypeB::LSecantBAlgorithm< Real >::dsrch(), ROL::LineSearch< Real >::getInitialAlpha(), ROL::LineSearch_U< Real >::getInitialAlpha(), ROL::CompositeStep< Real >::initialize(), ROL::InteriorPointStep< Real >::initialize(), ROL::InteriorPointStep< Real >::initialize(), ROL::PrimalDualActiveSetStep< Real >::initialize(), ROL::Step< Real >::initialize(), ROL::TrustRegionStep< Real >::initialize(), ROL::TypeB::ColemanLiAlgorithm< Real >::initialize(), ROL::TypeB::GradientAlgorithm< Real >::initialize(), ROL::TypeB::KelleySachsAlgorithm< Real >::initialize(), ROL::TypeB::LinMoreAlgorithm< Real >::initialize(), ROL::TypeB::LSecantBAlgorithm< Real >::initialize(), ROL::TypeB::NewtonKrylovAlgorithm< Real >::initialize(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::initialize(), ROL::TypeB::QuasiNewtonAlgorithm< Real >::initialize(), ROL::TypeB::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::initialize(), ROL::TypeP::InexactNewtonAlgorithm< Real >::initialize(), ROL::TypeP::iPianoAlgorithm< Real >::initialize(), ROL::TypeP::ProxGradientAlgorithm< Real >::initialize(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::initialize(), ROL::TypeP::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeP::TrustRegionAlgorithm< Real >::initialize(), ROL::TypeU::BundleAlgorithm< Real >::initialize(), ROL::TypeU::LineSearchAlgorithm< Real >::initialize(), ROL::TRUtils::initialRadius(), ROL::BackTracking< Real >::run(), ROL::BackTracking_U< Real >::run(), ROL::Bisection< Real >::run(), ROL::CubicInterp< Real >::run(), ROL::CubicInterp_U< Real >::run(), ROL::GoldenSection< Real >::run(), ROL::IterationScaling< Real >::run(), ROL::IterationScaling_U< Real >::run(), ROL::PathBasedTargetLevel< Real >::run(), ROL::PathBasedTargetLevel_U< Real >::run(), ROL::TypeB::ColemanLiAlgorithm< Real >::run(), ROL::TypeB::GradientAlgorithm< Real >::run(), ROL::TypeB::KelleySachsAlgorithm< Real >::run(), ROL::TypeB::LinMoreAlgorithm< Real >::run(), ROL::TypeB::LSecantBAlgorithm< Real >::run(), ROL::TypeB::NewtonKrylovAlgorithm< Real >::run(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::run(), ROL::TypeB::QuasiNewtonAlgorithm< Real >::run(), ROL::TypeB::SpectralGradientAlgorithm< Real >::run(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::run(), ROL::TypeP::InexactNewtonAlgorithm< Real >::run(), ROL::TypeP::iPianoAlgorithm< Real >::run(), ROL::TypeP::ProxGradientAlgorithm< Real >::run(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::run(), ROL::TypeP::SpectralGradientAlgorithm< Real >::run(), ROL::TypeP::TrustRegionAlgorithm< Real >::run(), ROL::TypeU::BundleAlgorithm< Real >::run(), ROL::TypeU::LineSearchAlgorithm< Real >::run(), ROL::LineSearch< Real >::status(), ROL::AugmentedLagrangianStep< Real >::update(), ROL::CompositeStep< Real >::update(), ROL::GradientStep< Real >::update(), ROL::NewtonKrylovStep< Real >::update(), ROL::NewtonStep< Real >::update(), ROL::NonlinearCGStep< Real >::update(), ROL::PrimalDualActiveSetStep< Real >::update(), ROL::ProjectedNewtonKrylovStep< Real >::update(), ROL::ProjectedNewtonStep< Real >::update(), ROL::ProjectedSecantStep< Real >::update(), ROL::SecantStep< Real >::update(), and ROL::TrustRegion< Real >::update().

◆ update() [2/2]

template<typename Real>
virtual void ROL::Objective< Real >::update ( const Vector< Real > & x,
bool flag = true,
int iter = -1 )
inlinevirtual

Update objective function.

This function updates the objective function at new iterations.

Parameters
[in]xis the new iterate.
[in]flagis true if the iterate has changed.
[in]iteris the outer algorithm iterations count.

Reimplemented in CLExactModel< Real >, and Objective_PoissonInversion< Real >.

Definition at line 75 of file ROL_Objective.hpp.

References ROL_UNUSED.

◆ value()

template<typename Real>
virtual Real ROL::Objective< Real >::value ( const Vector< Real > & x,
Real & tol )
pure virtual

Compute value.

This function returns the objective function value.

Parameters
[in]xis the current iterate.
[in]tolis a tolerance for inexact objective function computation.

Implemented in CLExactModel< Real >, CLTestObjective< Real >, NullObjective< Real >, Objective_BurgersControl< Real >, Objective_GrossPitaevskii< Real >, Objective_GrossPitaevskii< Real >, Objective_PoissonInversion< Real >, and Zakharov_Sacado_Objective< Real >.

Referenced by ROL::CompositeStep< Real >::accept(), ROL::BundleStep< Real >::compute(), ROL::CompositeStep< Real >::compute(), ROL::TypeB::LinMoreAlgorithm< Real >::computeValue(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::computeValue(), ROL::TypeP::TrustRegionAlgorithm< Real >::computeValue(), ROL::TypeP::TrustRegionAlgorithm< Real >::dbls(), ROL::TypeP::TrustRegionAlgorithm< Real >::dcauchy(), ROL::LineSearch_U< Real >::dirDeriv(), ROL::TypeP::TrustRegionAlgorithm< Real >::dncg(), ROL::TypeP::TrustRegionAlgorithm< Real >::dspg(), ROL::TypeP::TrustRegionAlgorithm< Real >::dspg2(), ROL::TypeB::LSecantBAlgorithm< Real >::dsrch(), ROL::LineSearch< Real >::getInitialAlpha(), ROL::LineSearch_U< Real >::getInitialAlpha(), ROL::CompositeStep< Real >::initialize(), ROL::InteriorPointStep< Real >::initialize(), ROL::InteriorPointStep< Real >::initialize(), ROL::PrimalDualActiveSetStep< Real >::initialize(), ROL::Step< Real >::initialize(), ROL::TrustRegionStep< Real >::initialize(), ROL::TypeB::ColemanLiAlgorithm< Real >::initialize(), ROL::TypeB::GradientAlgorithm< Real >::initialize(), ROL::TypeB::KelleySachsAlgorithm< Real >::initialize(), ROL::TypeB::LinMoreAlgorithm< Real >::initialize(), ROL::TypeB::LSecantBAlgorithm< Real >::initialize(), ROL::TypeB::NewtonKrylovAlgorithm< Real >::initialize(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::initialize(), ROL::TypeB::QuasiNewtonAlgorithm< Real >::initialize(), ROL::TypeB::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::initialize(), ROL::TypeP::InexactNewtonAlgorithm< Real >::initialize(), ROL::TypeP::iPianoAlgorithm< Real >::initialize(), ROL::TypeP::ProxGradientAlgorithm< Real >::initialize(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::initialize(), ROL::TypeP::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeP::TrustRegionAlgorithm< Real >::initialize(), ROL::TypeU::BundleAlgorithm< Real >::initialize(), ROL::TypeU::LineSearchAlgorithm< Real >::initialize(), ROL::TRUtils::initialRadius(), ROL::BackTracking< Real >::run(), ROL::BackTracking_U< Real >::run(), ROL::Bisection< Real >::run(), ROL::CubicInterp< Real >::run(), ROL::CubicInterp_U< Real >::run(), ROL::GoldenSection< Real >::run(), ROL::IterationScaling< Real >::run(), ROL::IterationScaling_U< Real >::run(), ROL::PathBasedTargetLevel< Real >::run(), ROL::PathBasedTargetLevel_U< Real >::run(), ROL::TypeB::ColemanLiAlgorithm< Real >::run(), ROL::TypeB::GradientAlgorithm< Real >::run(), ROL::TypeB::KelleySachsAlgorithm< Real >::run(), ROL::TypeB::NewtonKrylovAlgorithm< Real >::run(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::run(), ROL::TypeB::QuasiNewtonAlgorithm< Real >::run(), ROL::TypeB::SpectralGradientAlgorithm< Real >::run(), ROL::TypeP::InexactNewtonAlgorithm< Real >::run(), ROL::TypeP::iPianoAlgorithm< Real >::run(), ROL::TypeP::ProxGradientAlgorithm< Real >::run(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::run(), ROL::TypeP::SpectralGradientAlgorithm< Real >::run(), ROL::TypeU::BundleAlgorithm< Real >::run(), ROL::TypeU::LineSearchAlgorithm< Real >::run(), ROL::CompositeStep< Real >::update(), ROL::GradientStep< Real >::update(), ROL::NewtonKrylovStep< Real >::update(), ROL::NewtonStep< Real >::update(), ROL::NonlinearCGStep< Real >::update(), ROL::PrimalDualActiveSetStep< Real >::update(), ROL::ProjectedNewtonKrylovStep< Real >::update(), ROL::ProjectedNewtonStep< Real >::update(), ROL::ProjectedSecantStep< Real >::update(), ROL::SecantStep< Real >::update(), and ROL::TrustRegion< Real >::update().

◆ gradient()

template<typename Real>
virtual void ROL::Objective< Real >::gradient ( Vector< Real > & g,
const Vector< Real > & x,
Real & tol )
virtual

Compute gradient.

This function returns the objective function gradient.

Parameters
[out]gis the gradient.
[in]xis the current iterate.
[in]tolis a tolerance for inexact objective function computation.

The default implementation is a finite-difference approximation based on the function value. This requires the definition of a basis \(\{\phi_i\}\) for the optimization vectors x and the definition of a basis \(\{\psi_j\}\) for the dual optimization vectors (gradient vectors g). The bases must be related through the Riesz map, i.e., \( R \{\phi_i\} = \{\psi_j\}\), and this must be reflected in the implementation of the ROL::Vector::dual() method.

Reimplemented in CLExactModel< Real >, CLTestObjective< Real >, NullObjective< Real >, Objective_BurgersControl< Real >, Objective_GrossPitaevskii< Real >, Objective_GrossPitaevskii< Real >, Objective_PoissonInversion< Real >, and Zakharov_Sacado_Objective< Real >.

Referenced by ROL::CompositeStep< Real >::accept(), ROL::BundleStep< Real >::compute(), ROL::CompositeStep< Real >::compute(), ROL::PrimalDualActiveSetStep< Real >::computeCriticalityMeasure(), ROL::TypeB::LinMoreAlgorithm< Real >::computeGradient(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::computeGradient(), ROL::TypeP::TrustRegionAlgorithm< Real >::computeGradient(), ROL::CompositeStep< Real >::initialize(), ROL::InteriorPointStep< Real >::initialize(), ROL::InteriorPointStep< Real >::initialize(), ROL::Step< Real >::initialize(), ROL::TypeB::ColemanLiAlgorithm< Real >::initialize(), ROL::TypeB::GradientAlgorithm< Real >::initialize(), ROL::TypeB::KelleySachsAlgorithm< Real >::initialize(), ROL::TypeB::LSecantBAlgorithm< Real >::initialize(), ROL::TypeB::NewtonKrylovAlgorithm< Real >::initialize(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::initialize(), ROL::TypeB::QuasiNewtonAlgorithm< Real >::initialize(), ROL::TypeB::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeP::InexactNewtonAlgorithm< Real >::initialize(), ROL::TypeP::iPianoAlgorithm< Real >::initialize(), ROL::TypeP::ProxGradientAlgorithm< Real >::initialize(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::initialize(), ROL::TypeP::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeU::BundleAlgorithm< Real >::initialize(), ROL::TypeU::LineSearchAlgorithm< Real >::initialize(), ROL::TypeB::ColemanLiAlgorithm< Real >::run(), ROL::TypeB::GradientAlgorithm< Real >::run(), ROL::TypeB::KelleySachsAlgorithm< Real >::run(), ROL::TypeB::LSecantBAlgorithm< Real >::run(), ROL::TypeB::NewtonKrylovAlgorithm< Real >::run(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::run(), ROL::TypeB::QuasiNewtonAlgorithm< Real >::run(), ROL::TypeB::SpectralGradientAlgorithm< Real >::run(), ROL::TypeP::InexactNewtonAlgorithm< Real >::run(), ROL::TypeP::iPianoAlgorithm< Real >::run(), ROL::TypeP::ProxGradientAlgorithm< Real >::run(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::run(), ROL::TypeP::SpectralGradientAlgorithm< Real >::run(), ROL::TypeU::BundleAlgorithm< Real >::run(), ROL::TypeU::LineSearchAlgorithm< Real >::run(), ROL::LineSearch< Real >::status(), ROL::CompositeStep< Real >::update(), ROL::GradientStep< Real >::update(), ROL::NewtonKrylovStep< Real >::update(), ROL::NewtonStep< Real >::update(), ROL::NonlinearCGStep< Real >::update(), ROL::ProjectedNewtonKrylovStep< Real >::update(), ROL::ProjectedNewtonStep< Real >::update(), ROL::ProjectedSecantStep< Real >::update(), ROL::SecantStep< Real >::update(), ROL::TrustRegion< Real >::update(), and ROL::TrustRegionStep< Real >::updateGradient().

◆ dirDeriv()

template<typename Real>
virtual Real ROL::Objective< Real >::dirDeriv ( const Vector< Real > & x,
const Vector< Real > & d,
Real & tol )
virtual

Compute directional derivative.

This function returns the directional derivative of the objective function in the \(d\) direction.

Parameters
[in]xis the current iterate.
[in]dis the direction.
[in]tolis a tolerance for inexact objective function computation.

Referenced by ROL::LineSearch_U< Real >::dirDeriv(), and ROL::StdObjective< Real >::dirDeriv().

◆ hessVec()

template<typename Real>
virtual void ROL::Objective< Real >::hessVec ( Vector< Real > & hv,
const Vector< Real > & v,
const Vector< Real > & x,
Real & tol )
virtual

◆ invHessVec()

template<typename Real>
virtual void ROL::Objective< Real >::invHessVec ( Vector< Real > & hv,
const Vector< Real > & v,
const Vector< Real > & x,
Real & tol )
inlinevirtual

Apply inverse Hessian approximation to vector.

This function applies the inverse Hessian of the objective function to the vector \(v\).

Parameters
[out]hvis the action of the inverse Hessian on \(v\).
[in]vis the direction vector.
[in]xis the current iterate.
[in]tolis a tolerance for inexact objective function computation.

Definition at line 131 of file ROL_Objective.hpp.

References ROL_UNUSED.

Referenced by ROL::Newton_U< Real >::compute(), ROL::NewtonStep< Real >::compute(), ROL::ProjectedNewtonStep< Real >::compute(), ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionModel_U< Real >::validate().

◆ precond()

template<typename Real>
virtual void ROL::Objective< Real >::precond ( Vector< Real > & Pv,
const Vector< Real > & v,
const Vector< Real > & x,
Real & tol )
inlinevirtual

Apply preconditioner to vector.

This function applies a preconditioner for the Hessian of the objective function to the vector \(v\).

Parameters
[out]Pvis the action of the Hessian preconditioner on \(v\).
[in]vis the direction vector.
[in]xis the current iterate.
[in]tolis a tolerance for inexact objective function computation.

Definition at line 149 of file ROL_Objective.hpp.

References ROL::Vector< Real >::dual(), ROL_UNUSED, and ROL::Vector< Real >::set().

◆ prox()

template<typename Real>
virtual void ROL::Objective< Real >::prox ( Vector< Real > & Pv,
const Vector< Real > & v,
Real t,
Real & tol )
inlinevirtual

Compute the proximity operator.

This function returns the proximity operator.

Parameters
[out]Pvis the proximity operator applied to \(v\) (primal optimization vector).
[in]vis the input to the proximity operator (primal optimization vector).
[in]tis the proximity operator parameter (positive scalar).
[in]tolis a tolerance for inexact objective function computation.

Definition at line 163 of file ROL_Objective.hpp.

References ROL_UNUSED.

Referenced by ROL::TypeP::TrustRegionAlgorithm< Real >::dprox(), ROL::TypeP::InexactNewtonAlgorithm< Real >::initialize(), ROL::TypeP::iPianoAlgorithm< Real >::initialize(), ROL::TypeP::ProxGradientAlgorithm< Real >::initialize(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::initialize(), ROL::TypeP::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeP::TrustRegionAlgorithm< Real >::initialize(), ROL::TypeP::Algorithm< Real >::pgstep(), and ROL::TypeP::iPianoAlgorithm< Real >::run().

◆ proxJacVec()

template<typename Real>
virtual void ROL::Objective< Real >::proxJacVec ( Vector< Real > & Jv,
const Vector< Real > & v,
const Vector< Real > & x,
Real t,
Real & tol )
virtual

Apply the Jacobian of the proximity operator.

This function applies the Jacobian of the proximity operator.

Parameters
[out]Jvis the Jacobian of the proximity operator at \(x\) applied to \(v\) (primal optimization vector).
[in]vis the direction vector (primal optimization vector).
[in]xis the input to the proximity operator (primal optimization vector).
[in]tis the proximity operator parameter (positive scalar).
[in]tolis a tolerance for inexact objective function computation.

◆ checkGradient() [1/4]

template<typename Real>
virtual std::vector< std::vector< Real > > ROL::Objective< 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 )
inlinevirtual

Finite-difference gradient check.

This function computes a sequence of one-sided finite-difference checks for the gradient.
At each step of the sequence, the finite difference step size is decreased. The output compares the error

\[ \left| \frac{f(x+td) - f(x)}{t} - \langle \nabla f(x),d\rangle_{\mathcal{X}^*,\mathcal{X}}\right|. \]

if the approximation is first order. More generally, difference approximation is

\[ \frac{1}{t} \sum\limits_{i=1}^m w_i f(x+t c_i d) \]

where m = order+1, \(w_i\) are the difference weights and \(c_i\) are the difference steps

Parameters
[in]xis an optimization variable.
[in]dis a direction vector.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.
[in]numStepsis a parameter which dictates the number of finite difference steps.
[in]orderis the order of the finite difference approximation (1,2,3,4)

Definition at line 204 of file ROL_Objective.hpp.

References checkGradient(), ROL::Vector< Real >::dual(), and ROL_NUM_CHECKDERIV_STEPS.

Referenced by checkGradient(), checkGradient(), and main().

◆ checkGradient() [2/4]

template<typename Real>
virtual std::vector< std::vector< Real > > ROL::Objective< Real >::checkGradient ( const Vector< Real > & x,
const Vector< Real > & g,
const Vector< Real > & d,
const bool printToStream = true,
std::ostream & outStream = std::cout,
const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
const int order = 1 )
virtual

Finite-difference gradient check.

This function computes a sequence of one-sided finite-difference checks for the gradient.
At each step of the sequence, the finite difference step size is decreased. The output compares the error

\[ \left| \frac{f(x+td) - f(x)}{t} - \langle \nabla f(x),d\rangle_{\mathcal{X}^*,\mathcal{X}}\right|. \]

if the approximation is first order. More generally, difference approximation is

\[ \frac{1}{t} \sum\limits_{i=1}^m w_i f(x+t c_i d) \]

where m = order+1, \(w_i\) are the difference weights and \(c_i\) are the difference steps

Parameters
[in]xis an optimization variable.
[in]gis used to create a temporary gradient vector.
[in]dis a direction vector.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.
[in]numStepsis a parameter which dictates the number of finite difference steps.
[in]orderis the order of the finite difference approximation (1,2,3,4)

References ROL_NUM_CHECKDERIV_STEPS.

◆ checkGradient() [3/4]

template<typename Real>
virtual std::vector< std::vector< Real > > ROL::Objective< 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 )
inlinevirtual

Finite-difference gradient check with specified step sizes.

This function computes a sequence of one-sided finite-difference checks for the gradient.
At each step of the sequence, the finite difference step size is decreased. The output compares the error

\[ \left| \frac{f(x+td) - f(x)}{t} - \langle \nabla f(x),d\rangle_{\mathcal{X}^*,\mathcal{X}}\right|. \]

if the approximation is first order. More generally, difference approximation is

\[ \frac{1}{t} \sum\limits_{i=1}^m w_i f(x+t c_i d) \]

where m = order+1, \(w_i\) are the difference weights and \(c_i\) are the difference steps

Parameters
[in]xis an optimization variable.
[in]dis a direction vector.
[in]stepsis vector of steps of user-specified size.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.
[in]orderis the order of the finite difference approximation (1,2,3,4)

Definition at line 264 of file ROL_Objective.hpp.

References checkGradient(), and ROL::Vector< Real >::dual().

◆ checkGradient() [4/4]

template<typename Real>
virtual std::vector< std::vector< Real > > ROL::Objective< Real >::checkGradient ( const Vector< Real > & x,
const Vector< Real > & g,
const Vector< Real > & d,
const std::vector< Real > & steps,
const bool printToStream = true,
std::ostream & outStream = std::cout,
const int order = 1 )
virtual

Finite-difference gradient check with specified step sizes.

This function computes a sequence of one-sided finite-difference checks for the gradient.
At each step of the sequence, the finite difference step size is decreased. The output compares the error

\[ \left| \frac{f(x+td) - f(x)}{t} - \langle \nabla f(x),d\rangle_{\mathcal{X}^*,\mathcal{X}}\right|. \]

if the approximation is first order. More generally, difference approximation is

\[ \frac{1}{t} \sum\limits_{i=1}^m w_i f(x+t c_i d) \]

where m = order+1, \(w_i\) are the difference weights and \(c_i\) are the difference steps

Parameters
[in]xis an optimization variable.
[in]gis used to create a temporary gradient vector.
[in]dis a direction vector.
[in]stepsis vector of steps of user-specified size.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.
[in]orderis the order of the finite difference approximation (1,2,3,4)

◆ checkHessVec() [1/4]

template<typename Real>
virtual std::vector< std::vector< Real > > ROL::Objective< 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 )
inlinevirtual

Finite-difference Hessian-applied-to-vector check.

This function computes a sequence of one-sided finite-difference checks for the Hessian.
At each step of the sequence, the finite difference step size is decreased. The output compares the error

\[ \left\| \frac{\nabla f(x+tv) - \nabla f(x)}{t} - \nabla^2 f(x)v\right\|_{\mathcal{X}^*}, \]

if the approximation is first order. More generally, difference approximation is

\[ \frac{1}{t} \sum\limits_{i=1}^m w_i \nabla f(x+t c_i v), \]

where m = order+1, \(w_i\) are the difference weights and \(c_i\) are the difference steps

Parameters
[in]xis an optimization variable.
[in]vis a direction vector.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.
[in]numStepsis a parameter which dictates the number of finite difference steps.
[in]orderis the order of the finite difference approximation (1,2,3,4)

Definition at line 326 of file ROL_Objective.hpp.

References checkHessVec(), ROL::Vector< Real >::dual(), and ROL_NUM_CHECKDERIV_STEPS.

Referenced by checkHessVec(), checkHessVec(), and main().

◆ checkHessVec() [2/4]

template<typename Real>
virtual std::vector< std::vector< Real > > ROL::Objective< Real >::checkHessVec ( const Vector< Real > & x,
const Vector< Real > & hv,
const Vector< Real > & v,
const bool printToStream = true,
std::ostream & outStream = std::cout,
const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
const int order = 1 )
virtual

Finite-difference Hessian-applied-to-vector check.

This function computes a sequence of one-sided finite-difference checks for the Hessian.
At each step of the sequence, the finite difference step size is decreased. The output compares the error

\[ \left\| \frac{\nabla f(x+tv) - \nabla f(x)}{t} - \nabla^2 f(x)v\right\|_{\mathcal{X}^*}, \]

if the approximation is first order. More generally, difference approximation is

\[ \frac{1}{t} \sum\limits_{i=1}^m w_i \nabla f(x+t c_i v), \]

where m = order+1, \(w_i\) are the difference weights and \(c_i\) are the difference steps

Parameters
[in]xis an optimization variable.
[in]hvis used to create temporary gradient and Hessian-times-vector vectors.
[in]vis a direction vector.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.
[in]numStepsis a parameter which dictates the number of finite difference steps.
[in]orderis the order of the finite difference approximation (1,2,3,4)

References ROL_NUM_CHECKDERIV_STEPS.

◆ checkHessVec() [3/4]

template<typename Real>
virtual std::vector< std::vector< Real > > ROL::Objective< 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 )
inlinevirtual

Finite-difference Hessian-applied-to-vector check with specified step sizes.

This function computes a sequence of one-sided finite-difference checks for the Hessian.
At each step of the sequence, the finite difference step size is decreased. The output compares the error

\[ \left\| \frac{\nabla f(x+tv) - \nabla f(x)}{t} - \nabla^2 f(x)v\right\|_{\mathcal{X}^*}, \]

if the approximation is first order. More generally, difference approximation is

\[ \frac{1}{t} \sum\limits_{i=1}^m w_i \nabla f(x+t c_i v), \]

where m = order+1, \(w_i\) are the difference weights and \(c_i\) are the difference steps

Parameters
[in]xis an optimization variable.
[in]vis a direction vector.
[in]stepsis vector of steps of user-specified size.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.
[in]orderis the order of the finite difference approximation (1,2,3,4)

Definition at line 387 of file ROL_Objective.hpp.

References checkHessVec(), and ROL::Vector< Real >::dual().

◆ checkHessVec() [4/4]

template<typename Real>
virtual std::vector< std::vector< Real > > ROL::Objective< Real >::checkHessVec ( const Vector< Real > & x,
const Vector< Real > & hv,
const Vector< Real > & v,
const std::vector< Real > & steps,
const bool printToStream = true,
std::ostream & outStream = std::cout,
const int order = 1 )
virtual

Finite-difference Hessian-applied-to-vector check with specified step sizes.

This function computes a sequence of one-sided finite-difference checks for the Hessian.
At each step of the sequence, the finite difference step size is decreased. The output compares the error

\[ \left\| \frac{\nabla f(x+tv) - \nabla f(x)}{t} - \nabla^2 f(x)v\right\|_{\mathcal{X}^*}, \]

if the approximation is first order. More generally, difference approximation is

\[ \frac{1}{t} \sum\limits_{i=1}^m w_i \nabla f(x+t c_i v), \]

where m = order+1, \(w_i\) are the difference weights and \(c_i\) are the difference steps

Parameters
[in]xis an optimization variable.
[in]hvis used to create temporary gradient and Hessian-times-vector vectors.
[in]vis a direction vector.
[in]stepsis vector of steps of user-specified size.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.
[in]orderis the order of the finite difference approximation (1,2,3,4)

◆ checkHessSym() [1/2]

template<typename Real>
virtual std::vector< Real > ROL::Objective< Real >::checkHessSym ( const Vector< Real > & x,
const Vector< Real > & v,
const Vector< Real > & w,
const bool printToStream = true,
std::ostream & outStream = std::cout )
inlinevirtual

Hessian symmetry check.

This function checks the symmetry of the Hessian by comparing

\[ \langle \nabla^2f(x)v,w\rangle_{\mathcal{X}^*,\mathcal{X}} \quad\text{and}\quad \langle \nabla^2f(x)w,v\rangle_{\mathcal{X}^*,\mathcal{X}}. \]

Parameters
[in]xis an optimization variable.
[in]vis a direction vector.
[in]wis a direction vector.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.

Definition at line 442 of file ROL_Objective.hpp.

References checkHessSym(), and ROL::Vector< Real >::dual().

Referenced by checkHessSym().

◆ checkHessSym() [2/2]

template<typename Real>
virtual std::vector< Real > ROL::Objective< Real >::checkHessSym ( const Vector< Real > & x,
const Vector< Real > & hv,
const Vector< Real > & v,
const Vector< Real > & w,
const bool printToStream = true,
std::ostream & outStream = std::cout )
virtual

Hessian symmetry check.

This function checks the symmetry of the Hessian by comparing

\[ \langle \nabla^2f(x)v,w\rangle_{\mathcal{X}^*,\mathcal{X}} \quad\text{and}\quad \langle \nabla^2f(x)w,v\rangle_{\mathcal{X}^*,\mathcal{X}}. \]

Parameters
[in]xis an optimization variable.
[in]hvis used to create temporary Hessian-times-vector vectors.
[in]vis a direction vector.
[in]wis a direction vector.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.

◆ checkProxJacVec()

template<typename Real>
virtual std::vector< std::vector< Real > > ROL::Objective< 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 )
virtual

Finite-difference proximity operator Jacobian-applied-to-vector check.

This function computes a sequence of one-sided finite-difference checks for the proximity operator Jacobian.
At each step of the sequence, the finite difference step size is decreased. The output compares the error

\[ \left\| \frac{\mathrm{prox}_{t f}(x+tv) - \mathrm{prox}_{t f}(x)}{t} - J_{t f}(x+tv)v\right\|_{\mathcal{X}}, \]

if the approximation is first order. Note that in some cases the proximity operator is semismooth, which motivates the evaluation of \(J_{t f}\) at \(x+tv\).

Parameters
[in]xis an optimization vector.
[in]vis a direction vector.
[in]tis the proximity operator parameter.
[in]printToStreamis a flag that turns on/off output.
[out]outStreamis the output stream.
[in]numStepsis a parameter which dictates the number of finite difference steps.

References ROL_NUM_CHECKDERIV_STEPS.

◆ getParameter()

◆ setParameter()

Member Data Documentation

◆ prim_

template<typename Real>
Ptr<Vector<Real> > ROL::Objective< Real >::prim_
private

Definition at line 47 of file ROL_Objective.hpp.

Referenced by Objective().

◆ dual_

template<typename Real>
Ptr<Vector<Real> > ROL::Objective< Real >::dual_
private

Definition at line 47 of file ROL_Objective.hpp.

Referenced by Objective().

◆ basis_

template<typename Real>
Ptr<Vector<Real> > ROL::Objective< Real >::basis_
private

Definition at line 47 of file ROL_Objective.hpp.

Referenced by Objective().

◆ param_

template<typename Real>
std::vector<Real> ROL::Objective< Real >::param_
private

Definition at line 501 of file ROL_Objective.hpp.

Referenced by getParameter(), and setParameter().


The documentation for this class was generated from the following file: