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
63namespace ROL {
64
65template<typename Real>
66class TrustRegionModel_U : public Objective<Real> {
67private:
68 Ptr<Objective<Real>> obj_;
69 Ptr<const Vector<Real>> x_, g_;
70 Ptr<Vector<Real>> dual_;
71
72 Ptr<Secant<Real>> secant_;
75
76protected:
77 /***************************************************************************/
78 /********* BEGIN WRAPPERS FOR HESSIAN/PRECOND APPLICATION ****************/
79 /***************************************************************************/
80 void applyHessian(Vector<Real> &hv, const Vector<Real> &v, Real &tol) {
81 if ( useSecantHessVec_ && secant_ != nullPtr ) {
82 secant_->applyB(hv,v);
83 }
84 else {
85 obj_->hessVec(hv,v,*x_,tol);
86 }
87 }
88
89 void applyInvHessian(Vector<Real> &hv, const Vector<Real> &v, Real &tol) {
90 if ( useSecantHessVec_ && secant_ != nullPtr ) {
91 secant_->applyH(hv,v);
92 }
93 else {
94 obj_->invHessVec(hv,v,*x_,tol);
95 }
96 }
97
98 void applyPrecond(Vector<Real> &Pv, const Vector<Real> &v, Real &tol) {
99 if ( useSecantPrecond_ && secant_ != nullPtr ) {
100 secant_->applyH(Pv,v);
101 }
102 else {
103 obj_->precond(Pv,v,*x_,tol);
104 }
105 }
106 /***************************************************************************/
107 /********* END WRAPPERS FOR HESSIAN/PRECOND APPLICATION ******************/
108 /***************************************************************************/
109
110public:
111
113
114 TrustRegionModel_U(ParameterList &list,
115 const Ptr<Secant<Real>> &secant = nullPtr,
117 : obj_(nullPtr), x_(nullPtr), g_(nullPtr), secant_(secant) {
118 ParameterList &slist = list.sublist("General").sublist("Secant");
119 useSecantPrecond_ = slist.get("Use as Preconditioner", false);
120 useSecantHessVec_ = slist.get("Use as Hessian", false);
121 if (secant_ == nullPtr) {
122 secant_ = SecantFactory<Real>(list,mode);
123 }
124 }
125
126 void initialize(const Vector<Real> &x, const Vector<Real> &g) {
127 dual_ = g.clone();
128 }
129
130 // Some versions of Clang will issue a warning that update hides and
131 // overloaded virtual function without this using declaration
132 using Objective<Real>::update;
133
135 const Vector<Real> &x,
136 const Vector<Real> &g,
137 ETrustRegionU etr) {
138 if ( !useSecantHessVec_ &&
140 try {
141 Real htol = std::sqrt(ROL_EPSILON<Real>());
142 Ptr<Vector<Real>> v = g.clone();
143 Ptr<Vector<Real>> hv = x.clone();
144 obj.invHessVec(*hv,*v,x,htol);
145 }
146 catch (std::exception &e) {
147 useSecantHessVec_ = true;
148 }
149 }
150 }
151
152 virtual void setData(Objective<Real> &obj,
153 const Vector<Real> &x,
154 const Vector<Real> &g) {
155 obj_ = makePtrFromRef(obj);
156 x_ = makePtrFromRef(x);
157 g_ = makePtrFromRef(g);
158 }
159
160 void update(const Vector<Real> &x, const Vector<Real> &s,
161 const Vector<Real> &gold, const Vector<Real> &gnew,
162 const Real snorm, const int iter) {
163 // Update Secant Information
165 secant_->updateStorage(x,gnew,gold,s,snorm,iter);
166 }
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_->dot(s.dual());
177 return dual_->apply(s);
178 }
179
180 virtual void gradient( Vector<Real> &g, const Vector<Real> &s, Real &tol ) override {
181 applyHessian(g,s,tol);
182 g.plus(*g_);
183 }
184
185 virtual void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) override {
186 applyHessian(hv,v,tol);
187 }
188
189 virtual void invHessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) override {
190 applyInvHessian(hv,v,tol);
191 }
192
193 virtual void precond( Vector<Real> &Pv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) override {
194 applyPrecond(Pv,v,tol);
195 }
196 /***************************************************************************/
197 /********* END OBJECTIVE FUNCTION DEFINITIONS ****************************/
198 /***************************************************************************/
199
200 /***************************************************************************/
201 /********* BEGIN ACCESSOR FUNCTIONS **************************************/
202 /***************************************************************************/
203 virtual const Ptr<const Vector<Real>> getGradient(void) const {
204 return g_;
205 }
206
207 virtual const Ptr<const Vector<Real>> getIterate(void) const {
208 return x_;
209 }
210
211 virtual const Ptr<Objective<Real>> getObjective(void) const {
212 return obj_;
213 }
214 /***************************************************************************/
215 /********* END ACCESSOR FUNCTIONS ****************************************/
216 /***************************************************************************/
217
218}; // class TrustRegionModel_U
219
220} // namespace ROL
221
222
223#endif
Contains definitions of enums for trust region algorithms.
Provides the interface to evaluate objective functions.
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
Provides interface for and implements limited-memory secant operators.
Provides the interface to evaluate trust-region model functions.
void applyPrecond(Vector< Real > &Pv, const Vector< Real > &v, Real &tol)
virtual Real value(const Vector< Real > &s, Real &tol) override
Compute value.
void update(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &gold, const Vector< Real > &gnew, const Real snorm, const int iter)
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Apply Hessian approximation to vector.
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Apply inverse Hessian approximation to vector.
TrustRegionModel_U(ParameterList &list, const Ptr< Secant< Real > > &secant=nullPtr, ESecantMode mode=SECANTMODE_BOTH)
virtual const Ptr< const Vector< Real > > getIterate(void) const
virtual void gradient(Vector< Real > &g, const Vector< Real > &s, Real &tol) override
Compute gradient.
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Apply preconditioner to vector.
Ptr< const Vector< Real > > g_
virtual void setData(Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g)
Ptr< Objective< Real > > obj_
void validate(Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g, ETrustRegionU etr)
Ptr< const Vector< Real > > x_
void initialize(const Vector< Real > &x, const Vector< Real > &g)
virtual const Ptr< const Vector< Real > > getGradient(void) const
void applyHessian(Vector< Real > &hv, const Vector< Real > &v, Real &tol)
virtual const Ptr< Objective< Real > > getObjective(void) const
void applyInvHessian(Vector< Real > &hv, const Vector< Real > &v, Real &tol)
Defines the linear algebra or vector space interface.
virtual void plus(const Vector &x)=0
Compute , where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
ESecantMode
@ SECANTMODE_BOTH