ROL
Public Member Functions | Protected Member Functions | Private Attributes | List of all members
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. More...
 
virtual void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update objective function. More...
 
virtual Real value (const Vector< Real > &x, Real &tol)=0
 Compute value. More...
 
virtual void gradient (Vector< Real > &g, const Vector< Real > &x, Real &tol)
 Compute gradient. More...
 
virtual Real dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol)
 Compute directional derivative. More...
 
virtual void hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply Hessian approximation to vector. More...
 
virtual void invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply inverse Hessian approximation to vector. More...
 
virtual void precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply preconditioner to vector. More...
 
virtual void prox (Vector< Real > &Pv, const Vector< Real > &v, Real t, Real &tol)
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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:

It is strongly recommended that the user additionally overload:

The user may also overload:

Definition at line 78 of file ROL_Objective.hpp.

Constructor & Destructor Documentation

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

Definition at line 85 of file ROL_Objective.hpp.

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

Definition at line 87 of file ROL_Objective.hpp.

Member Function Documentation

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.

Reimplemented in ROL::MoreauYosidaObjective< Real >, ROL::ReducedDynamicObjective< Real >, ROL::InteriorPointObjective< Real >, ROL::Reduced_Objective_SimOpt< Real >, ROL::StochasticObjective< Real >, ROL::AugmentedLagrangianObjective< Real >, ROL::RiskNeutralObjective< Real >, ROL::ElasticObjective< Real >, ROL::FletcherObjectiveBase< Real >, ROL::ReducedDynamicStationaryControlsObjective< Real >, ROL::TypeP::InexactNewtonAlgorithm< Real >::NewtonObj, ROL::NonlinearLeastSquaresObjective< Real >, ROL::ChainRuleObjective< Real >, ROL::AffineTransformObjective< Real >, ROL::Objective_SimOpt< Real >, ROL::SimulatedObjectiveCVaR< Real >, ROL::CompositeObjective< Real >, ROL::ObjectiveFromConstraint< Real >, ROL::SlacklessObjective< Real >, ROL::StdObjective< Real >, ROL::SimulatedObjective< Real >, ROL::LinearCombinationObjective< Real >, ROL::ScaledObjective< Real >, ROL::MeanValueObjective< Real >, and ROL::RiskLessObjective< Real >.

Definition at line 96 of file ROL_Objective.hpp.

References ROL_UNUSED.

Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::accept(), ROL::CompositeStep< Real >::accept(), ROL::BundleStep< Real >::compute(), ROL::TypeU::TrustRegionAlgorithm< Real >::computeValue(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::computeValue(), ROL::TypeB::LinMoreAlgorithm< 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_U< Real >::getInitialAlpha(), ROL::LineSearch< Real >::getInitialAlpha(), ROL::TypeB::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeP::ProxGradientAlgorithm< Real >::initialize(), ROL::TypeP::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeB::GradientAlgorithm< Real >::initialize(), ROL::TypeP::iPianoAlgorithm< Real >::initialize(), ROL::TypeB::QuasiNewtonAlgorithm< Real >::initialize(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::initialize(), ROL::Step< Real >::initialize(), ROL::TypeU::BundleAlgorithm< Real >::initialize(), ROL::TypeB::LSecantBAlgorithm< Real >::initialize(), ROL::TypeU::TrustRegionAlgorithm< Real >::initialize(), ROL::TypeE::CompositeStepAlgorithm< Real >::initialize(), ROL::TypeB::ColemanLiAlgorithm< Real >::initialize(), ROL::TypeP::InexactNewtonAlgorithm< Real >::initialize(), ROL::TypeB::KelleySachsAlgorithm< Real >::initialize(), ROL::TypeU::LineSearchAlgorithm< Real >::initialize(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::initialize(), ROL::TypeB::LinMoreAlgorithm< Real >::initialize(), ROL::InteriorPointStep< Real >::initialize(), ROL::TypeB::NewtonKrylovAlgorithm< Real >::initialize(), ROL::TypeP::TrustRegionAlgorithm< Real >::initialize(), ROL::CompositeStep< Real >::initialize(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::initialize(), ROL::PrimalDualActiveSetStep< Real >::initialize(), ROL::TrustRegionStep< Real >::initialize(), ROL::TRUtils::initialRadius(), ROL::IterationScaling_U< Real >::run(), ROL::IterationScaling< Real >::run(), ROL::BackTracking_U< Real >::run(), ROL::CubicInterp_U< Real >::run(), ROL::CubicInterp< Real >::run(), ROL::BackTracking< Real >::run(), ROL::TypeB::SpectralGradientAlgorithm< Real >::run(), ROL::TypeB::GradientAlgorithm< Real >::run(), ROL::TypeP::iPianoAlgorithm< Real >::run(), ROL::TypeP::ProxGradientAlgorithm< Real >::run(), ROL::Bisection< Real >::run(), ROL::GoldenSection< Real >::run(), ROL::TypeP::SpectralGradientAlgorithm< Real >::run(), ROL::PathBasedTargetLevel_U< Real >::run(), ROL::PathBasedTargetLevel< Real >::run(), ROL::TypeB::QuasiNewtonAlgorithm< Real >::run(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::run(), ROL::TypeU::BundleAlgorithm< Real >::run(), ROL::TypeB::LSecantBAlgorithm< Real >::run(), ROL::TypeB::KelleySachsAlgorithm< Real >::run(), ROL::TypeU::TrustRegionAlgorithm< Real >::run(), ROL::TypeB::ColemanLiAlgorithm< Real >::run(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::run(), ROL::TypeB::LinMoreAlgorithm< Real >::run(), ROL::TypeU::LineSearchAlgorithm< Real >::run(), ROL::TypeP::InexactNewtonAlgorithm< Real >::run(), ROL::TypeB::NewtonKrylovAlgorithm< Real >::run(), ROL::TypeP::TrustRegionAlgorithm< Real >::run(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::run(), ROL::LineSearch< Real >::status(), ROL::NewtonStep< Real >::update(), ROL::GradientStep< Real >::update(), ROL::ProjectedNewtonStep< Real >::update(), ROL::NonlinearCGStep< Real >::update(), ROL::SecantStep< Real >::update(), ROL::ProjectedSecantStep< Real >::update(), ROL::TrustRegion< Real >::update(), ROL::NewtonKrylovStep< Real >::update(), ROL::ProjectedNewtonKrylovStep< Real >::update(), ROL::CompositeStep< Real >::update(), ROL::AugmentedLagrangianStep< Real >::update(), ROL::PrimalDualActiveSetStep< Real >::update(), and ROL::TypeE::CompositeStepAlgorithm< Real >::updateRadius().

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 ROL::BoundFletcher< Real >, Objective_PoissonInversion< Real >, ROL::ProjectedObjective< Real >, ROL::MoreauYosidaPenalty< Real >, ROL::ReducedDynamicObjective< Real >, ROL::Fletcher< Real >, ROL::InteriorPointPenalty< Real >, ROL::AugmentedLagrangian< Real >, ROL::StochasticObjective< Real >, ROL::Reduced_Objective_SimOpt< Real >, ROL::RiskNeutralObjective< Real >, CLExactModel< Real >, ROL::PH_bPOEObjective< Real >, ROL::PH_RiskObjective< Real >, ROL::PH_DeviationObjective< Real >, ROL::QuadraticPenalty< Real >, ROL::Reduced_AugmentedLagrangian_SimOpt< Real >, ROL::PH_ErrorObjective< Real >, ROL::PH_Objective< Real >, ROL::PH_RegretObjective< Real >, ROL::PH_ProbObjective< Real >, ROL::NonlinearLeastSquaresObjective_Dynamic< Real >, ROL::ReducedDynamicStationaryControlsObjective< Real >, ROL::ChainRuleObjective< Real >, ROL::InteriorPoint::PenalizedObjective< Real >, ROL::ObjectiveMMA< Real >, ROL::NonlinearLeastSquaresObjective< Real >, ROL::ConicApproximationModel< Real >, ROL::ConicApproximationModel< Real >, ROL::AffineTransformObjective< Real >, ROL::CompositeObjective< Real >, ROL::ObjectiveFromConstraint< Real >, ROL::SimulatedObjectiveCVaR< Real >, ROL::SlacklessObjective< Real >, ROL::Objective_SimOpt< Real >, ROL::LinearCombinationObjective< Real >, ROL::SimulatedObjective< Real >, ROL::StdObjective< Real >, ROL::MeanValueObjective< Real >, and ROL::RiskLessObjective< Real >.

Definition at line 109 of file ROL_Objective.hpp.

References ROL_UNUSED.

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 Objective_GrossPitaevskii< Real >, ROL::ZOO::Objective_PoissonInversion< Real >, ROL::ReducedDynamicObjective< Real >, ROL::BoundFletcher< Real >, Objective_BurgersControl< Real >, Objective_PoissonInversion< Real >, ROL::MoreauYosidaObjective< Real >, ROL::MoreauYosidaPenalty< Real >, ROL::ProjectedObjective< Real >, ROL::ColemanLiModel< Real >, ROL::InteriorPointObjective< Real >, ROL::Fletcher< Real >, ROL::InteriorPointPenalty< Real >, CLExactModel< Real >, ROL::AugmentedLagrangian< Real >, ROL::StochasticObjective< Real >, ROL::CDFObjective< Real >, ROL::RiskNeutralObjective< Real >, ROL::Reduced_Objective_SimOpt< Real >, ROL::ZOO::Objective_DiodeCircuit< Real >, ROL::AugmentedLagrangianObjective< Real >, ROL::TrustRegionModel_U< Real >, ROL::MomentObjective< Real >, ROL::PH_bPOEObjective< Real >, ROL::KelleySachsModel< Real >, ROL::PH_RiskObjective< Real >, Objective_GrossPitaevskii< Real >, ROL::FletcherObjectiveE< Real >, ROL::PH_DeviationObjective< Real >, ROL::TrustRegionModel< Real >, ROL::InteriorPoint::MeritFunction< Real >, ROL::ObjectiveFromBoundConstraint< Real >, ROL::QuadraticPenalty< Real >, ROL::ZOO::Objective_PoissonControl< Real >, Zakharov_Sacado_Objective< Real >, ROL::Reduced_AugmentedLagrangian_SimOpt< Real >, ROL::ObjectiveMMA< Real >, ROL::PH_ErrorObjective< Real >, ROL::PH_RegretObjective< Real >, ROL::PH_Objective< Real >, ROL::PH_ProbObjective< Real >, ROL::ChainRuleObjective< Real >, ROL::NonlinearLeastSquaresObjective_Dynamic< Real >, ROL::PointwiseCDFObjective< Real >, ROL::ReducedDynamicStationaryControlsObjective< Real >, ROL::ElasticObjective< Real >, ROL::InteriorPoint::PenalizedObjective< Real >, ROL::ZOO::Objective_Zakharov< Real >, ROL::ConicApproximationModel< Real >, ROL::ConicApproximationModel< Real >, ROL::TypeBIndicatorObjective< Real >, ROL::TypeP::InexactNewtonAlgorithm< Real >::NewtonObj, ROL::ZOO::Objective_HS45< Real >, ROL::NonlinearLeastSquaresObjective< Real >, ROL::ZOO::Objective_ParaboloidCircle< Real, XPrim, XDual >, ROL::ZOO::Objective_Rosenbrock< Real, XPrim, XDual >, ROL::ZOO::Objective_SimpleEqConstrained< Real, XPrim, XDual >, ROL::Objective_SimOpt< Real >, ROL::ZOO::Objective_HS2< Real >, ROL::ZOO::Objective_HS3< Real >, ROL::ZOO::Objective_HS38< Real >, ROL::ZOO::Objective_HS4< Real >, ROL::ZOO::Objective_HS5< Real >, ROL::AffineTransformObjective< Real >, ROL::ZOO::Objective_HS25< Real >, ROL::QuadraticObjective< Real >, ROL::l1Objective< Real >, ROL::ZOO::Objective_HS29< Real >, ROL::PQNObjective< Real >, ROL::SimulatedObjectiveCVaR< Real >, ROL::ZOO::Objective_HS32< Real >, ROL::ZOO::Objective_HS39< Real >, CLTestObjective< Real >, ROL::ZOO::Objective_Beale< Real >, ROL::ZOO::Objective_BVP< Real >, ROL::ZOO::Minimax3< Real >, ROL::CompositeObjective< Real >, ROL::LinearObjective< Real >, ROL::ObjectiveFromConstraint< Real >, ROL::StdObjective< Real >, ROL::SimulatedObjective< Real >, ROL::ZOO::Minimax1< Real >, ROL::ZOO::Minimax2< Real >, ROL::ZOO::Objective_HS24< Real >, ROL::SlacklessObjective< Real >, ROL::ZOO::Objective_HS1< Real >, ROL::BallIndicatorObjective< Real >, ROL::ZOO::Objective_LeastSquares< Real >, ROL::ZOO::Objective_Powell< Real >, ROL::LinearCombinationObjective< Real >, ROL::ZOO::Objective_FreudensteinRoth< Real >, ROL::ZOO::Objective_SumOfSquares< Real >, NullObjective< Real >, ROL::ScaledObjective< Real >, ROL::MeanValueObjective< Real >, ROL::RiskLessObjective< Real >, ROL::Objective_FSsolver< Real >, and ROL::LogBarrierObjective< Real >.

Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::accept(), ROL::CompositeStep< Real >::accept(), ROL::BundleStep< Real >::compute(), ROL::CompositeStep< Real >::compute(), ROL::TypeE::CompositeStepAlgorithm< Real >::computeTrial(), ROL::RandVarFunctional< Real >::computeValue(), ROL::TypeU::TrustRegionAlgorithm< Real >::computeValue(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::computeValue(), ROL::TypeB::LinMoreAlgorithm< 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_U< Real >::getInitialAlpha(), ROL::LineSearch< Real >::getInitialAlpha(), ROL::TypeB::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeP::ProxGradientAlgorithm< Real >::initialize(), ROL::TypeP::SpectralGradientAlgorithm< Real >::initialize(), ROL::TypeB::GradientAlgorithm< Real >::initialize(), ROL::TypeP::iPianoAlgorithm< Real >::initialize(), ROL::TypeB::QuasiNewtonAlgorithm< Real >::initialize(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::initialize(), ROL::Step< Real >::initialize(), ROL::TypeU::BundleAlgorithm< Real >::initialize(), ROL::TypeB::LSecantBAlgorithm< Real >::initialize(), ROL::TypeU::TrustRegionAlgorithm< Real >::initialize(), ROL::TypeE::CompositeStepAlgorithm< Real >::initialize(), ROL::TypeB::ColemanLiAlgorithm< Real >::initialize(), ROL::TypeP::InexactNewtonAlgorithm< Real >::initialize(), ROL::TypeB::KelleySachsAlgorithm< Real >::initialize(), ROL::TypeU::LineSearchAlgorithm< Real >::initialize(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::initialize(), ROL::TypeB::LinMoreAlgorithm< Real >::initialize(), ROL::InteriorPointStep< Real >::initialize(), ROL::TypeB::NewtonKrylovAlgorithm< Real >::initialize(), ROL::TypeP::TrustRegionAlgorithm< Real >::initialize(), ROL::CompositeStep< Real >::initialize(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::initialize(), ROL::PrimalDualActiveSetStep< Real >::initialize(), ROL::TrustRegionStep< Real >::initialize(), ROL::TRUtils::initialRadius(), ROL::IterationScaling_U< Real >::run(), ROL::IterationScaling< Real >::run(), ROL::CubicInterp_U< Real >::run(), ROL::BackTracking_U< Real >::run(), ROL::CubicInterp< Real >::run(), ROL::BackTracking< Real >::run(), ROL::TypeB::SpectralGradientAlgorithm< Real >::run(), ROL::TypeB::GradientAlgorithm< Real >::run(), ROL::TypeP::iPianoAlgorithm< Real >::run(), ROL::TypeP::ProxGradientAlgorithm< Real >::run(), ROL::TypeP::SpectralGradientAlgorithm< Real >::run(), ROL::Bisection< Real >::run(), ROL::GoldenSection< Real >::run(), ROL::PathBasedTargetLevel_U< Real >::run(), ROL::PathBasedTargetLevel< Real >::run(), ROL::TypeB::QuasiNewtonAlgorithm< Real >::run(), ROL::TypeP::QuasiNewtonAlgorithm< Real >::run(), ROL::TypeU::BundleAlgorithm< Real >::run(), ROL::TypeB::KelleySachsAlgorithm< Real >::run(), ROL::TypeB::ColemanLiAlgorithm< Real >::run(), ROL::TypeU::LineSearchAlgorithm< Real >::run(), ROL::TypeP::InexactNewtonAlgorithm< Real >::run(), ROL::TypeB::NewtonKrylovAlgorithm< Real >::run(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::run(), ROL::NewtonStep< Real >::update(), ROL::GradientStep< Real >::update(), ROL::ProjectedNewtonStep< Real >::update(), ROL::NonlinearCGStep< Real >::update(), ROL::SecantStep< Real >::update(), ROL::ProjectedSecantStep< Real >::update(), ROL::TrustRegion< Real >::update(), ROL::NewtonKrylovStep< Real >::update(), ROL::FletcherStep< Real >::update(), ROL::ProjectedNewtonKrylovStep< Real >::update(), ROL::CompositeStep< Real >::update(), ROL::PrimalDualActiveSetStep< Real >::update(), and ROL::TypeE::CompositeStepAlgorithm< Real >::updateRadius().

template<typename Real >
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 Objective_GrossPitaevskii< Real >, ROL::ReducedDynamicObjective< Real >, ROL::ZOO::Objective_PoissonInversion< Real >, ROL::BoundFletcher< Real >, Objective_BurgersControl< Real >, Objective_PoissonInversion< Real >, ROL::MoreauYosidaObjective< Real >, ROL::MoreauYosidaPenalty< Real >, ROL::InteriorPointPenalty< Real >, ROL::ProjectedObjective< Real >, ROL::InteriorPointObjective< Real >, ROL::Fletcher< Real >, ROL::ColemanLiModel< Real >, ROL::ObjectiveFromBoundConstraint< Real >, CLExactModel< Real >, ROL::AugmentedLagrangian< Real >, ROL::StochasticObjective< Real >, ROL::CDFObjective< Real >, ROL::RiskNeutralObjective< Real >, ROL::ZOO::Objective_DiodeCircuit< Real >, ROL::AugmentedLagrangianObjective< Real >, ROL::Reduced_Objective_SimOpt< Real >, ROL::ZOO::Objective_PoissonControl< Real >, Objective_GrossPitaevskii< Real >, ROL::MomentObjective< Real >, ROL::TrustRegionModel_U< Real >, ROL::KelleySachsModel< Real >, ROL::PH_bPOEObjective< Real >, ROL::PH_RiskObjective< Real >, ROL::PH_DeviationObjective< Real >, ROL::QuadraticPenalty< Real >, ROL::ObjectiveMMA< Real >, ROL::FletcherObjectiveE< Real >, ROL::TrustRegionModel< Real >, Zakharov_Sacado_Objective< Real >, ROL::Reduced_AugmentedLagrangian_SimOpt< Real >, ROL::Objective_SimOpt< Real >, ROL::PH_Objective< Real >, ROL::PH_ProbObjective< Real >, ROL::PH_ErrorObjective< Real >, ROL::PH_RegretObjective< Real >, ROL::ChainRuleObjective< Real >, ROL::PointwiseCDFObjective< Real >, ROL::ReducedDynamicStationaryControlsObjective< Real >, ROL::InteriorPoint::PenalizedObjective< Real >, ROL::NonlinearLeastSquaresObjective_Dynamic< Real >, ROL::ZOO::Objective_SimpleEqConstrained< Real, XPrim, XDual >, ROL::ZOO::Objective_Zakharov< Real >, ROL::ElasticObjective< Real >, ROL::SimulatedObjectiveCVaR< Real >, ROL::TypeP::InexactNewtonAlgorithm< Real >::NewtonObj, ROL::ZOO::Objective_Rosenbrock< Real, XPrim, XDual >, ROL::ZOO::Objective_ParaboloidCircle< Real, XPrim, XDual >, ROL::ZOO::Objective_HS45< Real >, ROL::ConicApproximationModel< Real >, ROL::ConicApproximationModel< Real >, ROL::ZOO::Objective_HS25< Real >, ROL::SimulatedObjective< Real >, ROL::ZOO::Objective_HS38< Real >, ROL::NonlinearLeastSquaresObjective< Real >, ROL::ZOO::Objective_HS2< Real >, ROL::ZOO::Objective_HS3< Real >, ROL::ZOO::Objective_HS4< Real >, ROL::ZOO::Objective_HS5< Real >, ROL::ZOO::Objective_LeastSquares< Real >, ROL::ZOO::Minimax3< Real >, ROL::ZOO::Objective_BVP< Real >, ROL::ZOO::Objective_HS29< Real >, ROL::l1Objective< Real >, ROL::ZOO::Objective_Beale< Real >, ROL::AffineTransformObjective< Real >, ROL::ZOO::Objective_HS32< Real >, ROL::QuadraticObjective< Real >, ROL::ZOO::Objective_HS39< Real >, CLTestObjective< Real >, ROL::ZOO::Minimax1< Real >, ROL::ZOO::Minimax2< Real >, ROL::PQNObjective< Real >, ROL::ZOO::Objective_HS24< Real >, ROL::ZOO::Objective_Powell< Real >, ROL::StdObjective< Real >, ROL::ZOO::Objective_FreudensteinRoth< Real >, ROL::ZOO::Objective_HS1< Real >, ROL::CompositeObjective< Real >, ROL::LinearObjective< Real >, ROL::ObjectiveFromConstraint< Real >, ROL::SlacklessObjective< Real >, ROL::LogBarrierObjective< Real >, ROL::ZOO::Objective_SumOfSquares< Real >, ROL::LinearCombinationObjective< Real >, ROL::ShiftedProxObjective< Real >, NullObjective< Real >, ROL::ScaledObjective< Real >, ROL::MeanValueObjective< Real >, ROL::RiskLessObjective< Real >, ROL::ZeroProxObjective< Real >, and ROL::Objective_FSsolver< Real >.

Definition at line 74 of file ROL_ObjectiveDef.hpp.

References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::basis(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::dimension(), ROL::Vector< Real >::dot(), ROL::Revert, ROL::Temp, ROL::update(), ROL::value, zero, and ROL::Vector< Real >::zero().

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

template<typename Real >
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.

Reimplemented in ROL::ProjectedObjective< Real >, ROL::ZOO::Objective_Zakharov< Real >, ROL::l1Objective< Real >, ROL::StdObjective< Real >, ROL::LogBarrierObjective< Real >, and ROL::SlacklessObjective< Real >.

Definition at line 54 of file ROL_ObjectiveDef.hpp.

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

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

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

Apply Hessian approximation to vector.

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

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

Reimplemented in Objective_GrossPitaevskii< Real >, ROL::ReducedDynamicObjective< Real >, ROL::BoundFletcher< Real >, Objective_BurgersControl< Real >, Objective_PoissonInversion< Real >, ROL::MoreauYosidaObjective< Real >, ROL::InteriorPointPenalty< Real >, ROL::MoreauYosidaPenalty< Real >, ROL::ObjectiveFromBoundConstraint< Real >, ROL::Fletcher< Real >, ROL::InteriorPointObjective< Real >, ROL::ProjectedObjective< Real >, ROL::ColemanLiModel< Real >, ROL::ZOO::Objective_DiodeCircuit< Real >, ROL::CDFObjective< Real >, ROL::Objective_SimOpt< Real >, ROL::RiskNeutralObjective< Real >, ROL::AugmentedLagrangian< Real >, CLExactModel< Real >, ROL::StochasticObjective< Real >, ROL::MomentObjective< Real >, ROL::AugmentedLagrangianObjective< Real >, Objective_GrossPitaevskii< Real >, ROL::Reduced_Objective_SimOpt< Real >, ROL::PH_bPOEObjective< Real >, ROL::KelleySachsModel< Real >, ROL::TrustRegionModel_U< Real >, ROL::PH_RiskObjective< Real >, ROL::ObjectiveMMA< Real >, ROL::QuadraticPenalty< Real >, ROL::PH_DeviationObjective< Real >, Zakharov_Sacado_Objective< Real >, ROL::TrustRegionModel< Real >, ROL::FletcherObjectiveE< Real >, ROL::Reduced_AugmentedLagrangian_SimOpt< Real >, ROL::PH_ProbObjective< Real >, ROL::InteriorPoint::PenalizedObjective< Real >, ROL::ChainRuleObjective< Real >, ROL::PH_Objective< Real >, ROL::ZOO::Objective_SimpleEqConstrained< Real, XPrim, XDual >, ROL::ReducedDynamicStationaryControlsObjective< Real >, ROL::PH_ErrorObjective< Real >, ROL::PH_RegretObjective< Real >, ROL::SimulatedObjective< Real >, ROL::ZOO::Objective_ParaboloidCircle< Real, XPrim, XDual >, ROL::NonlinearLeastSquaresObjective_Dynamic< Real >, ROL::TypeP::InexactNewtonAlgorithm< Real >::NewtonObj, ROL::ConicApproximationModel< Real >, ROL::ConicApproximationModel< Real >, ROL::ElasticObjective< Real >, ROL::ZOO::Objective_HS29< Real >, ROL::LogBarrierObjective< Real >, ROL::NonlinearLeastSquaresObjective< Real >, ROL::ZOO::Objective_HS32< Real >, ROL::ZOO::Objective_HS39< Real >, CLTestObjective< Real >, ROL::ZOO::Objective_HS24< Real >, ROL::StdObjective< Real >, ROL::AffineTransformObjective< Real >, ROL::QuadraticObjective< Real >, ROL::PQNObjective< Real >, ROL::CompositeObjective< Real >, ROL::LinearObjective< Real >, ROL::ObjectiveFromConstraint< Real >, ROL::SlacklessObjective< Real >, NullObjective< Real >, ROL::LinearCombinationObjective< Real >, ROL::ScaledObjective< Real >, ROL::MeanValueObjective< Real >, ROL::RiskLessObjective< Real >, and ROL::Objective_FSsolver< Real >.

Definition at line 95 of file ROL_ObjectiveDef.hpp.

References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::Revert, ROL::Vector< Real >::scale(), ROL::Temp, ROL::update(), zero, and ROL::Vector< Real >::zero().

Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::accept(), ROL::CompositeStep< Real >::accept(), ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::computeDenseHessian(), ROL::RandVarFunctional< Real >::computeHessVec(), ROL::computeScaledDenseHessian(), ROL::StdObjective< Real >::hessVec(), ROL::Reduced_Objective_SimOpt< Real >::hessVec(), ROL::ZOO::Objective_DiodeCircuit< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::hessVec(), ROL::TrustRegionStep< Real >::initialize(), main(), ROL::TypeB::PrimalDualActiveSetAlgorithm< Real >::run(), ROL::TypeE::CompositeStepAlgorithm< Real >::solveTangentialSubproblem(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().

template<typename Real>
virtual void ROL::Objective< Real >::invHessVec ( Vector< Real > &  hv,
const Vector< Real > &  v,
const Vector< Real > &  x,
Real &  tol 
)
inlinevirtual
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.

Reimplemented in ROL::ProjectedObjective< Real >, ROL::RiskNeutralObjective< Real >, ROL::StochasticObjective< Real >, ROL::KelleySachsModel< Real >, ROL::Reduced_Objective_SimOpt< Real >, ROL::TrustRegionModel_U< Real >, ROL::TrustRegionModel< Real >, ROL::ConicApproximationModel< Real >, ROL::ConicApproximationModel< Real >, ROL::NonlinearLeastSquaresObjective_Dynamic< Real >, ROL::StdObjective< Real >, ROL::NonlinearLeastSquaresObjective< Real >, ROL::SlacklessObjective< Real >, ROL::ScaledObjective< Real >, ROL::MeanValueObjective< Real >, and ROL::RiskLessObjective< Real >.

Definition at line 183 of file ROL_Objective.hpp.

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

template<typename Real>
virtual void ROL::Objective< Real >::prox ( Vector< Real > &  Pv,
const Vector< Real > &  v,
Real  t,
Real &  tol 
)
inlinevirtual
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 218 of file ROL_Objective.hpp.

References ROL::Vector< Real >::dual().

Referenced by ROL::Objective< Real >::checkGradient(), and main().

template<typename Real >
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)

Definition at line 119 of file ROL_ObjectiveDef.hpp.

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 278 of file ROL_Objective.hpp.

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

template<typename Real >
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)

Definition at line 138 of file ROL_ObjectiveDef.hpp.

References ROL::Vector< Real >::apply(), ROL::Vector< Real >::clone(), ROL::Finite_Difference_Arrays::shifts, ROL::Temp, ROL::update(), ROL::value, and ROL::Finite_Difference_Arrays::weights.

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 340 of file ROL_Objective.hpp.

References ROL::Vector< Real >::dual().

Referenced by ROL::Objective< Real >::checkHessVec(), and main().

template<typename Real >
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)

Definition at line 234 of file ROL_ObjectiveDef.hpp.

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 401 of file ROL_Objective.hpp.

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

template<typename Real >
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)

Definition at line 253 of file ROL_ObjectiveDef.hpp.

References ROL::Vector< Real >::clone(), ROL::Finite_Difference_Arrays::shifts, ROL::Temp, ROL::update(), and ROL::Finite_Difference_Arrays::weights.

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 456 of file ROL_Objective.hpp.

References ROL::Vector< Real >::dual().

Referenced by main().

template<typename Real >
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.

Definition at line 350 of file ROL_ObjectiveDef.hpp.

References ROL::Vector< Real >::apply(), ROL::Vector< Real >::clone(), ROL::Temp, and ROL::update().

template<typename Real>
const std::vector<Real> ROL::Objective< Real >::getParameter ( void  ) const
inlineprotected
template<typename Real>
virtual void ROL::Objective< Real >::setParameter ( const std::vector< Real > &  param)
inlinevirtual

Reimplemented in ROL::MoreauYosidaPenalty< Real >, ROL::CompositeObjective_SimOpt< Real >, ROL::Reduced_Objective_SimOpt< Real >, ROL::PH_bPOEObjective< Real >, ROL::PH_RiskObjective< Real >, ROL::PH_DeviationObjective< Real >, ROL::Reduced_AugmentedLagrangian_SimOpt< Real >, ROL::LinearCombinationObjective_SimOpt< Real >, ROL::PH_ProbObjective< Real >, ROL::PH_ErrorObjective< Real >, ROL::PH_RegretObjective< Real >, ROL::PH_Objective< Real >, ROL::NonlinearLeastSquaresObjective< Real >, ROL::AffineTransformObjective< Real >, ROL::CompositeObjective< Real >, ROL::SlacklessObjective< Real >, ROL::LinearCombinationObjective< Real >, ROL::RiskLessObjective< Real >, and ROL::ScaledObjective< Real >.

Definition at line 498 of file ROL_Objective.hpp.

References ROL::Objective< Real >::param_.

Referenced by ROL::RandVarFunctional< Real >::computeGradient(), ROL::RandVarFunctional< Real >::computeHessVec(), ROL::RandVarFunctional< Real >::computeValue(), ROL::ScaledObjective< Real >::setParameter(), ROL::RiskLessObjective< Real >::setParameter(), ROL::LinearCombinationObjective< Real >::setParameter(), ROL::SlacklessObjective< Real >::setParameter(), ROL::CompositeObjective< Real >::setParameter(), ROL::AffineTransformObjective< Real >::setParameter(), ROL::NonlinearLeastSquaresObjective< Real >::setParameter(), ROL::PH_Objective< Real >::setParameter(), ROL::PH_ErrorObjective< Real >::setParameter(), ROL::PH_RegretObjective< Real >::setParameter(), ROL::PH_ProbObjective< Real >::setParameter(), ROL::LinearCombinationObjective_SimOpt< Real >::setParameter(), ROL::Reduced_AugmentedLagrangian_SimOpt< Real >::setParameter(), ROL::PH_DeviationObjective< Real >::setParameter(), ROL::PH_RiskObjective< Real >::setParameter(), ROL::PH_bPOEObjective< Real >::setParameter(), ROL::Reduced_Objective_SimOpt< Real >::setParameter(), ROL::CompositeObjective_SimOpt< Real >::setParameter(), and ROL::MoreauYosidaPenalty< Real >::setParameter().

Member Data Documentation

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

Definition at line 81 of file ROL_Objective.hpp.

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

Definition at line 81 of file ROL_Objective.hpp.

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

Definition at line 81 of file ROL_Objective.hpp.

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

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