ROL
ROL_TrustRegionModel_U.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 #ifndef ROL_TRUSTREGIONMODEL_U_H
45 #define ROL_TRUSTREGIONMODEL_U_H
46 
47 #include "ROL_Objective.hpp"
48 #include "ROL_SecantFactory.hpp"
50 
63 namespace ROL {
64 
65 template<typename Real>
66 class TrustRegionModel_U : public Objective<Real> {
67 private:
68  Ptr<Objective<Real>> obj_;
69  Ptr<const Vector<Real>> x_, g_;
70  Ptr<Vector<Real>> dual_;
71  Real tol_;
72 
73  Ptr<Secant<Real>> secant_;
76 
77 protected:
78  /***************************************************************************/
79  /********* BEGIN WRAPPERS FOR HESSIAN/PRECOND APPLICATION ****************/
80  /***************************************************************************/
81  void applyHessian(Vector<Real> &hv, const Vector<Real> &v, Real &tol) {
82  if ( useSecantHessVec_ && secant_ != nullPtr ) {
83  secant_->applyB(hv,v);
84  }
85  else {
86  obj_->hessVec(hv,v,*x_,tol_);
87  }
88  }
89 
90  void applyInvHessian(Vector<Real> &hv, const Vector<Real> &v, Real &tol) {
91  if ( useSecantHessVec_ && secant_ != nullPtr ) {
92  secant_->applyH(hv,v);
93  }
94  else {
95  obj_->invHessVec(hv,v,*x_,tol_);
96  }
97  }
98 
99  void applyPrecond(Vector<Real> &Pv, const Vector<Real> &v, Real &tol) {
100  if ( useSecantPrecond_ && secant_ != nullPtr ) {
101  secant_->applyH(Pv,v);
102  }
103  else {
104  obj_->precond(Pv,v,*x_,tol_);
105  }
106  }
107  /***************************************************************************/
108  /********* END WRAPPERS FOR HESSIAN/PRECOND APPLICATION ******************/
109  /***************************************************************************/
110 
111 public:
112 
113  virtual ~TrustRegionModel_U() {}
114 
115  TrustRegionModel_U(ParameterList &list,
116  const Ptr<Secant<Real>> &secant = nullPtr,
118  : obj_(nullPtr), x_(nullPtr), g_(nullPtr), secant_(secant) {
119  ParameterList &slist = list.sublist("General").sublist("Secant");
120  useSecantPrecond_ = slist.get("Use as Preconditioner", false);
121  useSecantHessVec_ = slist.get("Use as Hessian", false);
122  if (secant_ == nullPtr) secant_ = SecantFactory<Real>(list,mode);
123  }
124 
125  void initialize(const Vector<Real> &x, const Vector<Real> &g) {
126  dual_ = g.clone();
127  }
128 
129  // Some versions of Clang will issue a warning that update hides and
130  // overloaded virtual function without this using declaration
132 
134  const Vector<Real> &x,
135  const Vector<Real> &g,
136  ETrustRegionU etr) {
137  if ( !useSecantHessVec_ &&
138  (etr == TRUSTREGION_U_DOGLEG || etr == TRUSTREGION_U_DOUBLEDOGLEG) ) {
139  try {
140  Real htol = std::sqrt(ROL_EPSILON<Real>());
141  Ptr<Vector<Real>> v = g.clone();
142  Ptr<Vector<Real>> hv = x.clone();
143  obj.invHessVec(*hv,*v,x,htol);
144  }
145  catch (std::exception &e) {
146  useSecantHessVec_ = true;
147  }
148  }
149  }
150 
151  virtual void setData(Objective<Real> &obj,
152  const Vector<Real> &x,
153  const Vector<Real> &g,
154  Real &tol) {
155  obj_ = makePtrFromRef(obj);
156  x_ = makePtrFromRef(x);
157  g_ = makePtrFromRef(g);
158  tol_ = tol;
159  }
160 
161  void update(const Vector<Real> &x, const Vector<Real> &s,
162  const Vector<Real> &gold, const Vector<Real> &gnew,
163  const Real snorm, const int iter) {
164  // Update Secant Information
166  secant_->updateStorage(x,gnew,gold,s,snorm,iter);
167  }
168 
169  /***************************************************************************/
170  /********* BEGIN OBJECTIVE FUNCTION DEFINITIONS **************************/
171  /***************************************************************************/
172  virtual Real value( const Vector<Real> &s, Real &tol ) override {
173  applyHessian(*dual_,s,tol);
174  dual_->scale(static_cast<Real>(0.5));
175  dual_->plus(*g_);
176  return dual_->apply(s);
177  }
178 
179  virtual void gradient( Vector<Real> &g, const Vector<Real> &s, Real &tol ) override {
180  applyHessian(g,s,tol);
181  g.plus(*g_);
182  }
183 
184  virtual void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) override {
185  applyHessian(hv,v,tol);
186  }
187 
188  virtual void invHessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) override {
189  applyInvHessian(hv,v,tol);
190  }
191 
192  virtual void precond( Vector<Real> &Pv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) override {
193  applyPrecond(Pv,v,tol);
194  }
195  /***************************************************************************/
196  /********* END OBJECTIVE FUNCTION DEFINITIONS ****************************/
197  /***************************************************************************/
198 
199  /***************************************************************************/
200  /********* BEGIN ACCESSOR FUNCTIONS **************************************/
201  /***************************************************************************/
202  virtual const Ptr<const Vector<Real>> getGradient(void) const {
203  return g_;
204  }
205 
206  virtual const Ptr<const Vector<Real>> getIterate(void) const {
207  return x_;
208  }
209 
210  virtual const Ptr<Objective<Real>> getObjective(void) const {
211  return obj_;
212  }
213  /***************************************************************************/
214  /********* END ACCESSOR FUNCTIONS ****************************************/
215  /***************************************************************************/
216 
217 }; // class TrustRegionModel_U
218 
219 } // namespace ROL
220 
221 
222 #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:80
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:57
Provides interface for and implements limited-memory secant operators.
Definition: ROL_Secant.hpp:79
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)