ROL
ROL_TrustRegionModel.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_H
45 #define ROL_TRUSTREGIONMODEL_H
46 
47 #include "ROL_Objective.hpp"
48 #include "ROL_BoundConstraint.hpp"
49 #include "ROL_Secant.hpp"
50 
64 namespace ROL {
65 
66 template <class Real>
67 class TrustRegionModel : public Objective<Real> {
68 private:
69  Ptr<Objective<Real>> obj_;
70  Ptr<BoundConstraint<Real>> bnd_;
71  Ptr<const Vector<Real>> x_, g_;
72  Ptr<Vector<Real>> dual_;
73  Ptr<Secant<Real>> secant_;
74 
75  const bool useSecantPrecond_;
76  const bool useSecantHessVec_;
77 
78  bool init_;
79 
80  void initialize(const Vector<Real> &s) {
81  if (!init_) {
82  dual_ = s.dual().clone();
83  init_ = true;
84  }
85  }
86 
87 protected:
88  /***************************************************************************/
89  /********* BEGIN WRAPPERS FOR HESSIAN/PRECOND APPLICATION ****************/
90  /***************************************************************************/
91  void applyHessian(Vector<Real> &hv, const Vector<Real> &v, Real &tol) {
92  if ( useSecantHessVec_ && secant_ != nullPtr ) {
93  secant_->applyB(hv,v);
94  }
95  else {
96  obj_->hessVec(hv,v,*x_,tol);
97  }
98  }
99 
100  void applyInvHessian(Vector<Real> &hv, const Vector<Real> &v, Real &tol) {
101  if ( useSecantHessVec_ && secant_ != nullPtr ) {
102  secant_->applyH(hv,v);
103  }
104  else {
105  obj_->invHessVec(hv,v,*x_,tol);
106  }
107  }
108 
109  void applyPrecond(Vector<Real> &Pv, const Vector<Real> &v, Real &tol) {
110  if ( useSecantPrecond_ && secant_ != nullPtr ) {
111  secant_->applyH(Pv,v);
112  }
113  else {
114  obj_->precond(Pv,v,*x_,tol);
115  }
116  }
117  /***************************************************************************/
118  /********* END WRAPPERS FOR HESSIAN/PRECOND APPLICATION ******************/
119  /***************************************************************************/
120 
121 public:
122 
123  virtual ~TrustRegionModel() {}
124 
126  const Vector<Real> &x, const Vector<Real> &g,
127  const Ptr<Secant<Real>> &secant = nullPtr,
128  const bool useSecantPrecond = false, const bool useSecantHessVec = false)
129  : obj_(makePtrFromRef(obj)), bnd_(makePtrFromRef(bnd)),
130  x_(makePtrFromRef(x)), g_(makePtrFromRef(g)),
131  secant_(secant), useSecantPrecond_(useSecantPrecond), useSecantHessVec_(useSecantHessVec),
132  init_(false) {}
133 
134  // Some versions of Clang will issue a warning that update hides and
135  // overloaded virtual function without this using declaration
137 
139  const Vector<Real> &x, const Vector<Real> &g,
140  const Ptr<Secant<Real>> &secant = nullPtr) {
141  obj_ = makePtrFromRef(obj);
142  bnd_ = makePtrFromRef(bnd);
143  x_ = makePtrFromRef(x);
144  g_ = makePtrFromRef(g);
145  secant_ = secant;
146  }
147 
148  /***************************************************************************/
149  /********* BEGIN OBJECTIVE FUNCTION DEFINITIONS **************************/
150  /***************************************************************************/
151  virtual Real value( const Vector<Real> &s, Real &tol ) {
152  initialize(s);
153  applyHessian(*dual_,s,tol);
154  dual_->scale(static_cast<Real>(0.5));
155  dual_->plus(*g_);
156  return dual_->dot(s.dual());
157  }
158 
159  virtual void gradient( Vector<Real> &g, const Vector<Real> &s, Real &tol ) {
160  applyHessian(g,s,tol);
161  g.plus(*g_);
162  }
163 
164  virtual void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) {
165  applyHessian(hv,v,tol);
166  }
167 
168  virtual void invHessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) {
169  applyInvHessian(hv,v,tol);
170  }
171 
172  virtual void precond( Vector<Real> &Pv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) {
173  applyPrecond(Pv,v,tol);
174  }
175  /***************************************************************************/
176  /********* END OBJECTIVE FUNCTION DEFINITIONS ****************************/
177  /***************************************************************************/
178 
179  /***************************************************************************/
180  /********* BEGIN ACCESSOR FUNCTIONS **************************************/
181  /***************************************************************************/
182  virtual const Ptr<const Vector<Real>> getGradient(void) const {
183  return g_;
184  }
185 
186  virtual const Ptr<const Vector<Real>> getIterate(void) const {
187  return x_;
188  }
189 
190  virtual const Ptr<Objective<Real>> getObjective(void) const {
191  return obj_;
192  }
193 
194  virtual const Ptr<BoundConstraint<Real>> getBoundConstraint(void) const {
195  if (!bnd_->isActivated()) {
196  return nullPtr;
197  }
198  return bnd_;
199  }
200  /***************************************************************************/
201  /********* END ACCESSOR FUNCTIONS ****************************************/
202  /***************************************************************************/
203 
204  virtual void dualTransform( Vector<Real> &tv, const Vector<Real> &v ) {
205  tv.set(v);
206  }
207 
208  virtual void primalTransform( Vector<Real> &tv, const Vector<Real> &v ) {
209  tv.set(v);
210  }
211 
212  virtual void updatePredictedReduction(Real &pred, const Vector<Real> &s) {}
213 
214  virtual void updateActualReduction(Real &ared, const Vector<Real> &s) {}
215 
216 }; // class TrustRegionModel
217 
218 } // namespace ROL
219 
220 
221 #endif
Provides the interface to evaluate objective functions.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: ROL_Vector.hpp:226
Ptr< Secant< Real > > secant_
virtual void plus(const Vector &x)=0
Compute , where .
Ptr< Objective< Real > > obj_
virtual void update(Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Ptr< Secant< Real >> &secant=nullPtr)
Ptr< const Vector< Real > > g_
Provides the interface to evaluate trust-region model functions.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
virtual const Ptr< const Vector< Real > > getGradient(void) const
virtual void updatePredictedReduction(Real &pred, const Vector< Real > &s)
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
void applyPrecond(Vector< Real > &Pv, const Vector< Real > &v, Real &tol)
void applyHessian(Vector< Real > &hv, const Vector< Real > &v, Real &tol)
TrustRegionModel(Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Ptr< Secant< Real >> &secant=nullPtr, const bool useSecantPrecond=false, const bool useSecantHessVec=false)
virtual const Ptr< Objective< Real > > getObjective(void) const
void applyInvHessian(Vector< Real > &hv, const Vector< Real > &v, Real &tol)
virtual Real value(const Vector< Real > &s, Real &tol)
Compute value.
Provides interface for and implements limited-memory secant operators.
Definition: ROL_Secant.hpp:70
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply Hessian approximation to vector.
void initialize(const Vector< Real > &s)
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply inverse Hessian approximation to vector.
virtual void gradient(Vector< Real > &g, const Vector< Real > &s, Real &tol)
Compute gradient.
Provides the interface to apply upper and lower bound constraints.
Ptr< BoundConstraint< Real > > bnd_
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
Ptr< Vector< Real > > dual_
virtual void updateActualReduction(Real &ared, const Vector< Real > &s)
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply preconditioner to vector.
Ptr< const Vector< Real > > x_
virtual const Ptr< BoundConstraint< Real > > getBoundConstraint(void) const
virtual void primalTransform(Vector< Real > &tv, const Vector< Real > &v)
virtual const Ptr< const Vector< Real > > getIterate(void) const