10 #ifndef ROL_DOUBLEDOGLEG_H
11 #define ROL_DOUBLEDOGLEG_H
26 ROL::Ptr<CauchyPoint<Real> >
cpt_;
28 ROL::Ptr<Vector<Real> >
s_;
29 ROL::Ptr<Vector<Real> >
v_;
30 ROL::Ptr<Vector<Real> >
Hp_;
38 cpt_ = ROL::makePtr<CauchyPoint<Real>>(parlist);
43 cpt_->initialize(x,s,g);
55 Real tol = std::sqrt(ROL_EPSILON<Real>());
56 const Real one(1),
zero(0), half(0.5), p2(0.2), p8(0.8), two(2);
62 Real sNnorm =
s_->norm();
63 Real tmp = -
s_->dot(s);
64 bool negCurv = (tmp >
zero ?
true :
false);
65 Real gsN = std::abs(tmp);
69 cpt_->run(s,snorm,iflag,iter,del,model);
86 Real gnorm = s.
norm();
87 Real gnorm2 = gnorm*gnorm;
89 Real gamma1 = gnorm/gBg;
90 Real gamma2 = gnorm/gsN;
91 Real eta = p8*gamma1*gamma2 + p2;
92 if (eta*sNnorm <= del || gBg <=
zero) {
101 if (gnorm2*gamma1 >= del) {
109 s.
scale(-gamma1*gnorm);
113 Real wNorm =
v_->dot(*
v_);
114 Real sigma = del*del-std::pow(gamma1*gnorm,two);
115 Real phi = s.
dot(*
v_);
116 Real theta = (-phi + std::sqrt(phi*phi+wNorm*sigma))/wNorm;
120 beta = (one-theta)*(-gamma1*gnorm);
124 pRed_ = -(alpha*(half*alpha-one)*gsN + half*beta*beta*gBg + beta*(one-alpha)*gnorm2);
ROL::Ptr< Vector< Real > > s_
void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual void scale(const Real alpha)=0
Compute where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
virtual void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
Contains definitions of custom data types in ROL.
Provides interface for and implements trust-region subproblem solvers.
Provides the interface to evaluate trust-region model functions.
Defines the linear algebra or vector space interface.
virtual const Ptr< const Vector< Real > > getGradient(void) const
virtual Real dot(const Vector &x) const =0
Compute where .
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
Provides interface for the double dog leg trust-region subproblem solver.
void setPredictedReduction(const Real pRed)
ROL::Ptr< Vector< Real > > v_
DoubleDogLeg(ROL::ParameterList &parlist)
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply Hessian approximation to vector.
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply inverse Hessian approximation to vector.
virtual void set(const Vector &x)
Set where .
virtual Real norm() const =0
Returns where .
void run(Vector< Real > &s, Real &snorm, int &iflag, int &iter, const Real del, TrustRegionModel< Real > &model)
ROL::Ptr< Vector< Real > > Hp_
ROL::Ptr< CauchyPoint< Real > > cpt_
virtual void primalTransform(Vector< Real > &tv, const Vector< Real > &v)