ROL
ROL_TrustRegionModel_U.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Rapid Optimization Library (ROL) Package
4 //
5 // Copyright 2014 NTESS and the ROL contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef ROL_TRUSTREGIONMODEL_U_H
11 #define ROL_TRUSTREGIONMODEL_U_H
12 
13 #include "ROL_Objective.hpp"
14 #include "ROL_SecantFactory.hpp"
16 
29 namespace ROL {
30 
31 template<typename Real>
32 class TrustRegionModel_U : public Objective<Real> {
33 private:
34  Ptr<Objective<Real>> obj_;
35  Ptr<const Vector<Real>> x_, g_;
36  Ptr<Vector<Real>> dual_;
37  Real tol_;
38 
39  Ptr<Secant<Real>> secant_;
42 
43 protected:
44  /***************************************************************************/
45  /********* BEGIN WRAPPERS FOR HESSIAN/PRECOND APPLICATION ****************/
46  /***************************************************************************/
47  void applyHessian(Vector<Real> &hv, const Vector<Real> &v, Real &tol) {
48  if ( useSecantHessVec_ && secant_ != nullPtr ) {
49  secant_->applyB(hv,v);
50  }
51  else {
52  obj_->hessVec(hv,v,*x_,tol_);
53  }
54  }
55 
56  void applyInvHessian(Vector<Real> &hv, const Vector<Real> &v, Real &tol) {
57  if ( useSecantHessVec_ && secant_ != nullPtr ) {
58  secant_->applyH(hv,v);
59  }
60  else {
61  obj_->invHessVec(hv,v,*x_,tol_);
62  }
63  }
64 
65  void applyPrecond(Vector<Real> &Pv, const Vector<Real> &v, Real &tol) {
66  if ( useSecantPrecond_ && secant_ != nullPtr ) {
67  secant_->applyH(Pv,v);
68  }
69  else {
70  obj_->precond(Pv,v,*x_,tol_);
71  }
72  }
73  /***************************************************************************/
74  /********* END WRAPPERS FOR HESSIAN/PRECOND APPLICATION ******************/
75  /***************************************************************************/
76 
77 public:
78 
79  virtual ~TrustRegionModel_U() {}
80 
81  TrustRegionModel_U(ParameterList &list,
82  const Ptr<Secant<Real>> &secant = nullPtr,
84  : obj_(nullPtr), x_(nullPtr), g_(nullPtr), secant_(secant) {
85  ParameterList &slist = list.sublist("General").sublist("Secant");
86  useSecantPrecond_ = slist.get("Use as Preconditioner", false);
87  useSecantHessVec_ = slist.get("Use as Hessian", false);
88  if (secant_ == nullPtr) secant_ = SecantFactory<Real>(list,mode);
89  }
90 
91  void initialize(const Vector<Real> &x, const Vector<Real> &g) {
92  dual_ = g.clone();
93  }
94 
95  // Some versions of Clang will issue a warning that update hides and
96  // overloaded virtual function without this using declaration
98 
100  const Vector<Real> &x,
101  const Vector<Real> &g,
102  ETrustRegionU etr) {
103  if ( !useSecantHessVec_ &&
104  (etr == TRUSTREGION_U_DOGLEG || etr == TRUSTREGION_U_DOUBLEDOGLEG) ) {
105  try {
106  Real htol = std::sqrt(ROL_EPSILON<Real>());
107  Ptr<Vector<Real>> v = g.clone();
108  Ptr<Vector<Real>> hv = x.clone();
109  obj.invHessVec(*hv,*v,x,htol);
110  }
111  catch (std::exception &e) {
112  useSecantHessVec_ = true;
113  }
114  }
115  }
116 
117  virtual void setData(Objective<Real> &obj,
118  const Vector<Real> &x,
119  const Vector<Real> &g,
120  Real &tol) {
121  obj_ = makePtrFromRef(obj);
122  x_ = makePtrFromRef(x);
123  g_ = makePtrFromRef(g);
124  tol_ = tol;
125  }
126 
127  void update(const Vector<Real> &x, const Vector<Real> &s,
128  const Vector<Real> &gold, const Vector<Real> &gnew,
129  const Real snorm, const int iter) {
130  // Update Secant Information
132  secant_->updateStorage(x,gnew,gold,s,snorm,iter);
133  }
134 
135  /***************************************************************************/
136  /********* BEGIN OBJECTIVE FUNCTION DEFINITIONS **************************/
137  /***************************************************************************/
138  virtual Real value( const Vector<Real> &s, Real &tol ) override {
139  applyHessian(*dual_,s,tol);
140  dual_->scale(static_cast<Real>(0.5));
141  dual_->plus(*g_);
142  return dual_->apply(s);
143  }
144 
145  virtual void gradient( Vector<Real> &g, const Vector<Real> &s, Real &tol ) override {
146  applyHessian(g,s,tol);
147  g.plus(*g_);
148  }
149 
150  virtual void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) override {
151  applyHessian(hv,v,tol);
152  }
153 
154  virtual void invHessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) override {
155  applyInvHessian(hv,v,tol);
156  }
157 
158  virtual void precond( Vector<Real> &Pv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) override {
159  applyPrecond(Pv,v,tol);
160  }
161  /***************************************************************************/
162  /********* END OBJECTIVE FUNCTION DEFINITIONS ****************************/
163  /***************************************************************************/
164 
165  /***************************************************************************/
166  /********* BEGIN ACCESSOR FUNCTIONS **************************************/
167  /***************************************************************************/
168  virtual const Ptr<const Vector<Real>> getGradient(void) const {
169  return g_;
170  }
171 
172  virtual const Ptr<const Vector<Real>> getIterate(void) const {
173  return x_;
174  }
175 
176  virtual const Ptr<Objective<Real>> getObjective(void) const {
177  return obj_;
178  }
179  /***************************************************************************/
180  /********* END ACCESSOR FUNCTIONS ****************************************/
181  /***************************************************************************/
182 
183 }; // class TrustRegionModel_U
184 
185 } // namespace ROL
186 
187 
188 #endif
Provides the interface to evaluate objective functions.
Ptr< Objective< Real > > obj_
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual Real value(const Vector< Real > &s, Real &tol) override
Compute value.
virtual void plus(const Vector &x)=0
Compute , where .
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Apply inverse Hessian approximation to vector.
Ptr< Secant< Real > > secant_
virtual void gradient(Vector< Real > &g, const Vector< Real > &s, Real &tol) override
Compute gradient.
Ptr< const Vector< Real > > x_
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:46
void validate(Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g, ETrustRegionU etr)
virtual const Ptr< const Vector< Real > > getIterate(void) const
void initialize(const Vector< Real > &x, const Vector< Real > &g)
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Apply Hessian approximation to vector.
Provides the interface to evaluate trust-region model functions.
Contains definitions of enums for trust region algorithms.
virtual const Ptr< const Vector< Real > > getGradient(void) const
ESecantMode
Definition: ROL_Secant.hpp:23
Provides interface for and implements limited-memory secant operators.
Definition: ROL_Secant.hpp:45
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
Ptr< const Vector< Real > > g_
virtual const Ptr< Objective< Real > > getObjective(void) const
void applyPrecond(Vector< Real > &Pv, const Vector< Real > &v, Real &tol)
void update(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &gold, const Vector< Real > &gnew, const Real snorm, const int iter)
ETrustRegionU
Enumeration of trust-region solver types.
void applyHessian(Vector< Real > &hv, const Vector< Real > &v, Real &tol)
TrustRegionModel_U(ParameterList &list, const Ptr< Secant< Real >> &secant=nullPtr, ESecantMode mode=SECANTMODE_BOTH)
void applyInvHessian(Vector< Real > &hv, const Vector< Real > &v, Real &tol)
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Apply preconditioner to vector.
virtual void setData(Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g, Real &tol)