44 #ifndef ROL_TRUNCATEDCG_H 
   45 #define ROL_TRUNCATEDCG_H 
   61   ROL::Ptr<Vector<Real> > 
s_;
 
   62   ROL::Ptr<Vector<Real> > 
g_;
 
   63   ROL::Ptr<Vector<Real> > 
v_;
 
   64   ROL::Ptr<Vector<Real> > 
p_;
 
   65   ROL::Ptr<Vector<Real> > 
Hp_;
 
   78     Real em4(1e-4), em2(1e-2);
 
   79     maxit_ = parlist.sublist(
"General").sublist(
"Krylov").get(
"Iteration Limit",20);
 
   80     tol1_  = parlist.sublist(
"General").sublist(
"Krylov").get(
"Absolute Tolerance",em4);
 
   81     tol2_  = parlist.sublist(
"General").sublist(
"Krylov").get(
"Relative Tolerance",em2);
 
  102     Real tol = std::sqrt(ROL_EPSILON<Real>());
 
  103     const Real 
zero(0), one(1), two(2), half(0.5);
 
  107     Real snorm2(0), s1norm2(0);
 
  110     Real gnorm = 
g_->norm(), normg = gnorm;
 
  111     const Real gtol = std::min(
tol1_,
tol2_*gnorm);
 
  115     p_->set(*
v_); 
p_->scale(-one);
 
  116     Real pnorm2 = 
v_->dot(
g_->dual());
 
  117     if ( pnorm2 <= 
zero ) {
 
  124     Real kappa(0), beta(0), sigma(0), alpha(0), tmp(0), sMp(0);
 
  125     Real gv = 
v_->dot(
g_->dual());
 
  128     for (iter = 0; iter < 
maxit_; iter++) {
 
  132       kappa = 
p_->dot(
Hp_->dual());
 
  134         sigma = (-sMp+sqrt(sMp*sMp+pnorm2*(del*del-snorm2)))/pnorm2;
 
  143       s1norm2 = snorm2 + two*alpha*sMp + alpha*alpha*pnorm2;
 
  145       if (s1norm2 >= del*del) {
 
  146         sigma = (-sMp+sqrt(sMp*sMp+pnorm2*(del*del-snorm2)))/pnorm2;
 
  152       pRed_ += half*alpha*gv;
 
  157       g_->axpy(alpha,*
Hp_);
 
  165       gv    = 
v_->dot(
g_->dual());
 
  170       sMp    = beta*(sMp+alpha*pnorm2);
 
  171       pnorm2 = gv + beta*beta*pnorm2; 
 
  175       pRed_ += sigma*(gv-half*sigma*kappa);
 
  178     if (iter == maxit_) {
 
  194     Real tol = std::sqrt(ROL_EPSILON<Real>());
 
  196     const Real gtol = std::min(
tol1_,
tol2_*gnorm);
 
  200     ROL::Ptr<Vector<Real> > sc = x.
clone();
 
  201     cauchypoint(*sc,scnorm,del,iflag,iter,x,grad,gnorm,pObj);
 
  202     ROL::Ptr<Vector<Real> > xc = x.
clone();
 
  209     Real snorm2  = snorm*snorm;
 
  210     ROL::Ptr<Vector<Real> > s1 = x.
clone();
 
  215     ROL::Ptr<Vector<Real> > g = x.
clone(); 
 
  217     ROL::Ptr<Vector<Real> > Hs = x.
clone();
 
  220     Real normg = g->norm();
 
  223     ROL::Ptr<Vector<Real> > v  = x.
clone();
 
  227     ROL::Ptr<Vector<Real> > p = x.
clone(); 
 
  230     Real pnorm2 = v->dot(*g);
 
  233     ROL::Ptr<Vector<Real> > Hp = x.
clone();
 
  242     Real gv     = v->dot(*g);
 
  246     for (iter = 0; iter < 
maxit_; iter++) {
 
  251         sigma = (-sMp+sqrt(sMp*sMp+pnorm2*(del*del-snorm2)))/pnorm2;
 
  260       s1norm2 = snorm2 + 2.0*alpha*sMp + alpha*alpha*pnorm2;
 
  262       if (s1norm2 >= del*del) {
 
  263         sigma = (-sMp+sqrt(sMp*sMp+pnorm2*(del*del-snorm2)))/pnorm2;
 
  269       pRed_ += 0.5*alpha*gv;
 
  287       sMp    = beta*(sMp+alpha*pnorm2);
 
  288       pnorm2 = gv + beta*beta*pnorm2; 
 
  291       pRed_ += sigma*(gv-0.5*sigma*kappa);
 
  294     if (iter == maxit_) {
 
ROL::Ptr< Vector< Real > > v_
 
void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
 
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)
 
void reducedHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &d, const Vector< Real > &x, Real &tol)
Apply the reduced Hessian to a vector, v. The reduced Hessian first removes elements of v correspondi...
 
Contains definitions of custom data types in ROL. 
 
ROL::Ptr< Vector< Real > > s_
 
Provides interface for and implements trust-region subproblem solvers. 
 
Provides the interface to evaluate trust-region model functions. 
 
virtual void zero()
Set to zero vector. 
 
Defines the linear algebra or vector space interface. 
 
virtual const Ptr< const Vector< Real > > getGradient(void) const 
 
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
 
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
 
void setPredictedReduction(const Real pRed)
 
TruncatedCG(ROL::ParameterList &parlist)
 
ROL::Ptr< Vector< Real > > g_
 
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply Hessian approximation to vector. 
 
void run(Vector< Real > &s, Real &snorm, int &iflag, int &iter, const Real del, TrustRegionModel< Real > &model)
 
ROL::Ptr< Vector< Real > > Hp_
 
Provides interface for truncated CG trust-region subproblem solver. 
 
ROL::Ptr< Vector< Real > > primalVector_
 
virtual void set(const Vector &x)
Set  where . 
 
virtual Real norm() const =0
Returns  where . 
 
void reducedPrecond(Vector< Real > &Mv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &d, const Vector< Real > &x, Real &tol)
Apply the reduced preconditioner to a vector, v. The reduced preconditioner first removes elements of...
 
ROL::Ptr< Vector< Real > > p_
 
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply preconditioner to vector. 
 
virtual void primalTransform(Vector< Real > &tv, const Vector< Real > &v)