ROL
|
Provides the interface to evaluate interior trust-region model functions from the Coleman-Li bound constrained trust-region algorithm. More...
#include <ROL_ColemanLiModel.hpp>
Public Member Functions | |
ColemanLiModel (Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Real stepBackMax=0.9999, const Real stepBackScale=1.0, const bool singleReflect=true, const Ptr< Secant< Real >> &secant=nullPtr, const bool useSecantPrecond=false, const bool useSecantHessVec=false) | |
void | update (Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Ptr< Secant< Real >> &secant=nullPtr) |
void | setRadius (const Real del) |
Real | value (const Vector< Real > &s, Real &tol) |
Compute value. More... | |
void | gradient (Vector< Real > &g, const Vector< Real > &s, Real &tol) |
Compute gradient. More... | |
void | hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
Apply Hessian approximation to vector. More... | |
void | dualTransform (Vector< Real > &tv, const Vector< Real > &v) |
void | primalTransform (Vector< Real > &tiv, const Vector< Real > &v) |
void | updatePredictedReduction (Real &pred, const Vector< Real > &s) |
void | updateActualReduction (Real &ared, const Vector< Real > &s) |
Public Member Functions inherited from ROL::TrustRegionModel< Real > | |
virtual | ~TrustRegionModel () |
TrustRegionModel (Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Ptr< Secant< Real >> &secant=nullPtr, const bool useSecantPrecond=false, const bool useSecantHessVec=false) | |
virtual void | invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
Apply inverse Hessian approximation to vector. More... | |
virtual void | precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
Apply preconditioner to vector. More... | |
virtual const Ptr< const Vector< Real > > | getGradient (void) const |
virtual const Ptr< const Vector< Real > > | getIterate (void) const |
virtual const Ptr< Objective < Real > > | getObjective (void) const |
virtual const Ptr < BoundConstraint< Real > > | getBoundConstraint (void) const |
Public Member Functions inherited from ROL::Objective< Real > | |
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 | dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol) |
Compute directional derivative. 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 > ¶m) |
Private Member Functions | |
void | applyD (Vector< Real > &Dv, const Vector< Real > &v) |
void | applyInverseD (Vector< Real > &Dv, const Vector< Real > &v) |
void | applyC (Vector< Real > &Cv, const Vector< Real > &v) |
void | constructC (void) |
void | constructInverseD (void) |
void | initialize (const Vector< Real > &x, const Vector< Real > &g) |
void | getScalarBounds (Real &lowerBound, Real &upperBound, const Vector< Real > &p) |
Real | minimize1D (Real &tau, const Real lowerBound, const Real upperBound, const Vector< Real > &p) |
Real | computeCauchyPoint (void) |
void | computeReflectiveStep (Vector< Real > &Rv, const Vector< Real > &v, const Vector< Real > &Dv) |
void | computeFullReflectiveStep (Vector< Real > &Rv, const Vector< Real > &v, const Vector< Real > &Dv) |
Real | computeAlpha (const Vector< Real > &d) |
bool | isStrictlyFeasibleStep (const Vector< Real > &d) const |
Private Attributes | |
Ptr< Vector< Real > > | prim_ |
Ptr< Vector< Real > > | dual_ |
Ptr< Vector< Real > > | hv_ |
Ptr< Vector< Real > > | step_ |
Ptr< Vector< Real > > | cauchyStep_ |
Ptr< Vector< Real > > | cauchyScal_ |
Ptr< Vector< Real > > | reflectStep_ |
Ptr< Vector< Real > > | reflectScal_ |
Ptr< Vector< Real > > | Dmat_ |
Ptr< Vector< Real > > | Cmat_ |
Ptr< Vector< Real > > | lx_ |
Ptr< Vector< Real > > | ux_ |
Real | TRradius_ |
const Real | stepBackMax_ |
const Real | stepBackScale_ |
const bool | singleReflect_ |
Real | sCs_ |
Real | pred_ |
Elementwise::Multiply< Real > | mult_ |
Elementwise::Divide< Real > | div_ |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::TrustRegionModel< Real > | |
void | applyHessian (Vector< Real > &hv, const Vector< Real > &v, Real &tol) |
void | applyInvHessian (Vector< Real > &hv, const Vector< Real > &v, Real &tol) |
void | applyPrecond (Vector< Real > &Pv, const Vector< Real > &v, Real &tol) |
Protected Member Functions inherited from ROL::Objective< Real > | |
const std::vector< Real > | getParameter (void) const |
Provides the interface to evaluate interior trust-region model functions from the Coleman-Li bound constrained trust-region algorithm.
Definition at line 26 of file ROL_ColemanLiModel.hpp.
|
inline |
Definition at line 189 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::initialize().
|
inlineprivate |
Definition at line 45 of file ROL_ColemanLiModel.hpp.
References ROL::Vector< Real >::applyBinary(), ROL::ColemanLiModel< Real >::div_, ROL::ColemanLiModel< Real >::Dmat_, and ROL::Vector< Real >::set().
|
inlineprivate |
Definition at line 51 of file ROL_ColemanLiModel.hpp.
References ROL::Vector< Real >::applyBinary(), ROL::ColemanLiModel< Real >::Dmat_, ROL::ColemanLiModel< Real >::mult_, and ROL::Vector< Real >::set().
Referenced by ROL::ColemanLiModel< Real >::computeCauchyPoint(), ROL::ColemanLiModel< Real >::dualTransform(), ROL::ColemanLiModel< Real >::gradient(), ROL::ColemanLiModel< Real >::hessVec(), ROL::ColemanLiModel< Real >::minimize1D(), ROL::ColemanLiModel< Real >::primalTransform(), and ROL::ColemanLiModel< Real >::value().
|
inlineprivate |
Definition at line 57 of file ROL_ColemanLiModel.hpp.
References ROL::Vector< Real >::applyBinary(), ROL::ColemanLiModel< Real >::Cmat_, ROL::ColemanLiModel< Real >::mult_, and ROL::Vector< Real >::set().
Referenced by ROL::ColemanLiModel< Real >::hessVec(), and ROL::ColemanLiModel< Real >::primalTransform().
|
inlineprivate |
Definition at line 62 of file ROL_ColemanLiModel.hpp.
References ROL::apply, ROL::ColemanLiModel< Real >::Cmat_, ROL::TrustRegionModel< Real >::getBoundConstraint(), ROL::TrustRegionModel< Real >::getGradient(), ROL::ColemanLiModel< Real >::mult_, ROL::ColemanLiModel< Real >::prim_, and zero.
Referenced by ROL::ColemanLiModel< Real >::update().
|
inlineprivate |
Definition at line 98 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::Dmat_, ROL::TrustRegionModel< Real >::getBoundConstraint(), ROL::TrustRegionModel< Real >::getGradient(), ROL::TrustRegionModel< Real >::getIterate(), ROL::ColemanLiModel< Real >::mult_, ROL::ColemanLiModel< Real >::prim_, ROL::ColemanLiModel< Real >::reflectScal_, ROL::ColemanLiModel< Real >::reflectStep_, and zero.
Referenced by ROL::ColemanLiModel< Real >::update().
|
inlineprivate |
Definition at line 171 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::cauchyScal_, ROL::ColemanLiModel< Real >::cauchyStep_, ROL::Vector< Real >::clone(), ROL::ColemanLiModel< Real >::Cmat_, ROL::ColemanLiModel< Real >::Dmat_, ROL::ColemanLiModel< Real >::dual_, ROL::ColemanLiModel< Real >::hv_, ROL::ColemanLiModel< Real >::lx_, ROL::ColemanLiModel< Real >::prim_, ROL::ColemanLiModel< Real >::reflectScal_, ROL::ColemanLiModel< Real >::reflectStep_, ROL::ColemanLiModel< Real >::step_, and ROL::ColemanLiModel< Real >::ux_.
Referenced by ROL::ColemanLiModel< Real >::ColemanLiModel().
|
inlinevirtual |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 202 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::constructC(), ROL::ColemanLiModel< Real >::constructInverseD(), and ROL::TrustRegionModel< Real >::update().
|
inline |
Definition at line 210 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::TRradius_.
|
inlinevirtual |
Compute value.
This function returns the objective function value.
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 218 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::applyInverseD(), ROL::Vector< Real >::dual(), ROL::TrustRegionModel< Real >::getGradient(), ROL::ColemanLiModel< Real >::hessVec(), ROL::ColemanLiModel< Real >::hv_, and ROL::ColemanLiModel< Real >::prim_.
Referenced by ROL::ColemanLiModel< Real >::computeCauchyPoint(), and ROL::ColemanLiModel< Real >::primalTransform().
|
inlinevirtual |
Compute gradient.
This function returns the objective function gradient.
[out] | g | is the gradient. |
[in] | x | is the current iterate. |
[in] | tol | is 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 from ROL::TrustRegionModel< Real >.
Definition at line 230 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::applyInverseD(), ROL::TrustRegionModel< Real >::getGradient(), ROL::ColemanLiModel< Real >::hessVec(), ROL::Vector< Real >::plus(), and ROL::ColemanLiModel< Real >::prim_.
|
inlinevirtual |
Apply Hessian approximation to vector.
This function applies the Hessian of the objective function to the vector \(v\).
[out] | hv | is the the action of the Hessian on \(v\). |
[in] | v | is the direction vector. |
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 237 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::applyC(), ROL::TrustRegionModel< Real >::applyHessian(), ROL::ColemanLiModel< Real >::applyInverseD(), ROL::ColemanLiModel< Real >::dual_, ROL::TrustRegionModel< Real >::getGradient(), ROL::Vector< Real >::plus(), and ROL::ColemanLiModel< Real >::prim_.
Referenced by ROL::ColemanLiModel< Real >::gradient(), ROL::ColemanLiModel< Real >::minimize1D(), and ROL::ColemanLiModel< Real >::value().
|
inlinevirtual |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 251 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::applyInverseD().
|
inlinevirtual |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 255 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::applyC(), ROL::ColemanLiModel< Real >::applyInverseD(), ROL::ColemanLiModel< Real >::cauchyScal_, ROL::ColemanLiModel< Real >::cauchyStep_, ROL::ColemanLiModel< Real >::computeCauchyPoint(), ROL::ColemanLiModel< Real >::computeFullReflectiveStep(), ROL::ColemanLiModel< Real >::computeReflectiveStep(), ROL::ColemanLiModel< Real >::getScalarBounds(), ROL::ColemanLiModel< Real >::isStrictlyFeasibleStep(), ROL::ColemanLiModel< Real >::minimize1D(), ROL::ColemanLiModel< Real >::pred_, ROL::ColemanLiModel< Real >::prim_, ROL::ColemanLiModel< Real >::reflectScal_, ROL::ColemanLiModel< Real >::reflectStep_, ROL::Vector< Real >::scale(), ROL::ColemanLiModel< Real >::sCs_, ROL::Vector< Real >::set(), ROL::ColemanLiModel< Real >::singleReflect_, ROL::ColemanLiModel< Real >::step_, ROL::ColemanLiModel< Real >::stepBackMax_, ROL::ColemanLiModel< Real >::stepBackScale_, and ROL::ColemanLiModel< Real >::value().
|
inlinevirtual |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 341 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::pred_.
|
inlinevirtual |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 345 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::sCs_.
|
inlineprivate |
Definition at line 351 of file ROL_ColemanLiModel.hpp.
References ROL::apply, ROL::TrustRegionModel< Real >::getBoundConstraint(), ROL::TrustRegionModel< Real >::getIterate(), ROL::Vector< Real >::norm(), ROL::ColemanLiModel< Real >::prim_, and ROL::ColemanLiModel< Real >::TRradius_.
Referenced by ROL::ColemanLiModel< Real >::computeCauchyPoint(), and ROL::ColemanLiModel< Real >::primalTransform().
|
inlineprivate |
Definition at line 405 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::applyInverseD(), ROL::Vector< Real >::dual(), ROL::TrustRegionModel< Real >::getGradient(), ROL::ColemanLiModel< Real >::hessVec(), ROL::ColemanLiModel< Real >::hv_, and ROL::ColemanLiModel< Real >::prim_.
Referenced by ROL::ColemanLiModel< Real >::computeCauchyPoint(), and ROL::ColemanLiModel< Real >::primalTransform().
|
inlineprivate |
Definition at line 428 of file ROL_ColemanLiModel.hpp.
References ROL::ColemanLiModel< Real >::applyInverseD(), ROL::ColemanLiModel< Real >::cauchyScal_, ROL::ColemanLiModel< Real >::cauchyStep_, ROL::TrustRegionModel< Real >::getGradient(), ROL::ColemanLiModel< Real >::getScalarBounds(), ROL::ColemanLiModel< Real >::minimize1D(), and ROL::ColemanLiModel< Real >::value().
Referenced by ROL::ColemanLiModel< Real >::primalTransform().
|
inlineprivate |
Definition at line 451 of file ROL_ColemanLiModel.hpp.
References ROL::apply, ROL::Vector< Real >::applyBinary(), ROL::ColemanLiModel< Real >::computeAlpha(), ROL::TrustRegionModel< Real >::getIterate(), ROL::ColemanLiModel< Real >::mult_, ROL::ColemanLiModel< Real >::prim_, and ROL::Vector< Real >::set().
Referenced by ROL::ColemanLiModel< Real >::primalTransform().
|
inlineprivate |
Definition at line 477 of file ROL_ColemanLiModel.hpp.
References ROL::apply, ROL::Vector< Real >::applyBinary(), ROL::TrustRegionModel< Real >::getIterate(), ROL::ColemanLiModel< Real >::mult_, ROL::ColemanLiModel< Real >::prim_, and ROL::Vector< Real >::set().
Referenced by ROL::ColemanLiModel< Real >::primalTransform().
|
inlineprivate |
Definition at line 502 of file ROL_ColemanLiModel.hpp.
References ROL::apply, ROL::TrustRegionModel< Real >::getIterate(), ROL::ColemanLiModel< Real >::lx_, ROL::ColemanLiModel< Real >::ux_, and zero.
Referenced by ROL::ColemanLiModel< Real >::computeReflectiveStep().
|
inlineprivate |
Definition at line 535 of file ROL_ColemanLiModel.hpp.
References ROL::apply, ROL::TrustRegionModel< Real >::getIterate(), and ROL::ColemanLiModel< Real >::prim_.
Referenced by ROL::ColemanLiModel< Real >::primalTransform().
|
private |
Definition at line 28 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::computeFullReflectiveStep(), ROL::ColemanLiModel< Real >::computeReflectiveStep(), ROL::ColemanLiModel< Real >::constructC(), ROL::ColemanLiModel< Real >::constructInverseD(), ROL::ColemanLiModel< Real >::getScalarBounds(), ROL::ColemanLiModel< Real >::gradient(), ROL::ColemanLiModel< Real >::hessVec(), ROL::ColemanLiModel< Real >::initialize(), ROL::ColemanLiModel< Real >::isStrictlyFeasibleStep(), ROL::ColemanLiModel< Real >::minimize1D(), ROL::ColemanLiModel< Real >::primalTransform(), and ROL::ColemanLiModel< Real >::value().
|
private |
Definition at line 28 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::hessVec(), and ROL::ColemanLiModel< Real >::initialize().
|
private |
Definition at line 28 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::initialize(), ROL::ColemanLiModel< Real >::minimize1D(), and ROL::ColemanLiModel< Real >::value().
|
private |
Definition at line 29 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::initialize(), and ROL::ColemanLiModel< Real >::primalTransform().
|
private |
Definition at line 30 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::computeCauchyPoint(), ROL::ColemanLiModel< Real >::initialize(), and ROL::ColemanLiModel< Real >::primalTransform().
|
private |
Definition at line 30 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::computeCauchyPoint(), ROL::ColemanLiModel< Real >::initialize(), and ROL::ColemanLiModel< Real >::primalTransform().
|
private |
Definition at line 31 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::constructInverseD(), ROL::ColemanLiModel< Real >::initialize(), and ROL::ColemanLiModel< Real >::primalTransform().
|
private |
Definition at line 31 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::constructInverseD(), ROL::ColemanLiModel< Real >::initialize(), and ROL::ColemanLiModel< Real >::primalTransform().
|
private |
Definition at line 32 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::applyD(), ROL::ColemanLiModel< Real >::applyInverseD(), ROL::ColemanLiModel< Real >::constructInverseD(), and ROL::ColemanLiModel< Real >::initialize().
|
private |
Definition at line 33 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::applyC(), ROL::ColemanLiModel< Real >::constructC(), and ROL::ColemanLiModel< Real >::initialize().
|
private |
Definition at line 34 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::computeAlpha(), and ROL::ColemanLiModel< Real >::initialize().
|
private |
Definition at line 34 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::computeAlpha(), and ROL::ColemanLiModel< Real >::initialize().
|
private |
Definition at line 36 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::getScalarBounds(), and ROL::ColemanLiModel< Real >::setRadius().
|
private |
Definition at line 37 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::primalTransform().
|
private |
Definition at line 37 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::primalTransform().
|
private |
Definition at line 38 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::primalTransform().
|
private |
Definition at line 39 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::primalTransform(), and ROL::ColemanLiModel< Real >::updateActualReduction().
|
private |
Definition at line 39 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::primalTransform(), and ROL::ColemanLiModel< Real >::updatePredictedReduction().
|
private |
Definition at line 41 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::applyC(), ROL::ColemanLiModel< Real >::applyInverseD(), ROL::ColemanLiModel< Real >::computeFullReflectiveStep(), ROL::ColemanLiModel< Real >::computeReflectiveStep(), ROL::ColemanLiModel< Real >::constructC(), and ROL::ColemanLiModel< Real >::constructInverseD().
|
private |
Definition at line 42 of file ROL_ColemanLiModel.hpp.
Referenced by ROL::ColemanLiModel< Real >::applyD().