44#ifndef ROL_TRUSTREGIONMODEL_U_H
45#define ROL_TRUSTREGIONMODEL_U_H
65template<
typename Real>
69 Ptr<const Vector<Real>>
x_,
g_;
80 void applyHessian(Vector<Real> &hv,
const Vector<Real> &v, Real &tol) {
81 if ( useSecantHessVec_ && secant_ != nullPtr ) {
82 secant_->applyB(hv,v);
85 obj_->hessVec(hv,v,*x_,tol);
90 if ( useSecantHessVec_ && secant_ != nullPtr ) {
91 secant_->applyH(hv,v);
94 obj_->invHessVec(hv,v,*x_,tol);
98 void applyPrecond(Vector<Real> &Pv,
const Vector<Real> &v, Real &tol) {
99 if ( useSecantPrecond_ && secant_ != nullPtr ) {
100 secant_->applyH(Pv,v);
103 obj_->precond(Pv,v,*x_,tol);
116 ESecantMode
mode = SECANTMODE_BOTH)
117 :
obj_(nullPtr), x_(nullPtr), g_(nullPtr), secant_(secant) {
118 ParameterList &slist = list.sublist(
"General").sublist(
"Secant");
119 useSecantPrecond_ = slist.get(
"Use as Preconditioner",
false);
120 useSecantHessVec_ = slist.get(
"Use as Hessian",
false);
121 if (secant_ == nullPtr) {
122 secant_ = SecantFactory<Real>(list,
mode);
126 void initialize(
const Vector<Real> &x,
const Vector<Real> &g) {
132 using Objective<Real>::update;
135 const Vector<Real> &x,
136 const Vector<Real> &g,
138 if ( !useSecantHessVec_ &&
139 (etr == TRUSTREGION_U_DOGLEG || etr == TRUSTREGION_U_DOUBLEDOGLEG) ) {
141 Real htol = std::sqrt(ROL_EPSILON<Real>());
142 Ptr<Vector<Real>> v = g.clone();
143 Ptr<Vector<Real>> hv = x.clone();
144 obj.invHessVec(*hv,*v,x,htol);
146 catch (std::exception &e) {
147 useSecantHessVec_ =
true;
153 const Vector<Real> &x,
154 const Vector<Real> &g) {
155 obj_ = makePtrFromRef(obj);
156 x_ = makePtrFromRef(x);
157 g_ = makePtrFromRef(g);
160 void update(
const Vector<Real> &x,
const Vector<Real> &s,
161 const Vector<Real> &gold,
const Vector<Real> &gnew,
164 if (useSecantHessVec_ || useSecantPrecond_) {
165 secant_->updateStorage(x,gnew,gold,s,
snorm,
iter);
172 virtual Real
value(
const Vector<Real> &s, Real &tol )
override {
174 dual_->scale(
static_cast<Real
>(0.5));
177 return dual_->apply(s);
180 virtual void gradient( Vector<Real> &g,
const Vector<Real> &s, Real &tol )
override {
185 virtual void hessVec( Vector<Real> &hv,
const Vector<Real> &v,
const Vector<Real> &s, Real &tol )
override {
189 virtual void invHessVec( Vector<Real> &hv,
const Vector<Real> &v,
const Vector<Real> &s, Real &tol )
override {
193 virtual void precond( Vector<Real> &Pv,
const Vector<Real> &v,
const Vector<Real> &s, Real &tol )
override {
207 virtual const Ptr<const Vector<Real>>
getIterate(
void)
const {
Contains definitions of enums for trust region algorithms.
Provides interface for and implements limited-memory secant operators.
Provides the interface to evaluate trust-region model functions.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Ptr< Objective< Real > > obj_
virtual const Ptr< const Vector< Real > > getIterate(void) const
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
void update(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &gold, const Vector< Real > &gnew, const Real snorm, const int iter)
void initialize(const Vector< Real > &x, const Vector< Real > &g)
Ptr< Vector< Real > > dual_
virtual ~TrustRegionModel_U()
virtual void setData(Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g)
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Ptr< const Vector< Real > > g_
Ptr< const Vector< Real > > x_
TrustRegionModel_U(ParameterList &list, const Ptr< Secant< Real > > &secant=nullPtr, ESecantMode mode=SECANTMODE_BOTH)
void applyPrecond(Vector< Real > &Pv, const Vector< Real > &v, Real &tol)
void applyHessian(Vector< Real > &hv, const Vector< Real > &v, Real &tol)
virtual const Ptr< Objective< Real > > getObjective(void) const
Ptr< Secant< Real > > secant_
virtual Real value(const Vector< Real > &s, Real &tol) override
virtual void gradient(Vector< Real > &g, const Vector< Real > &s, Real &tol) override
void applyInvHessian(Vector< Real > &hv, const Vector< Real > &v, Real &tol)
virtual const Ptr< const Vector< Real > > getGradient(void) const
void validate(Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g, ETrustRegionU etr)