ROL
ROL_HelperFunctions.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 
49 #ifndef ROL_HELPERFUNCTIONS_HPP
50 #define ROL_HELPERFUNCTIONS_HPP
51 
52 #include "ROL_Vector.hpp"
53 #include "ROL_Objective.hpp"
54 #include "ROL_BoundConstraint.hpp"
55 #include "ROL_Secant.hpp"
56 #include "Teuchos_SerialDenseMatrix.hpp"
57 #include "Teuchos_SerialDenseVector.hpp"
58 #include "Teuchos_LAPACK.hpp"
59 
60 namespace ROL {
61 
62  template<class Real>
63  Teuchos::SerialDenseMatrix<int, Real> computeDenseHessian(Objective<Real> &obj, const Vector<Real> &x) {
64 
65  Real tol = std::sqrt(ROL_EPSILON);
66 
67  int dim = x.dimension();
68  Teuchos::SerialDenseMatrix<int, Real> H(dim, dim);
69 
70  Teuchos::RCP<Vector<Real> > e = x.clone();
71  Teuchos::RCP<Vector<Real> > h = x.clone();
72 
73  for (int i=0; i<dim; i++) {
74  e = x.basis(i);
75  obj.hessVec(*h, *e, x, tol);
76  for (int j=0; j<dim; j++) {
77  e = x.basis(j);
78  H(j,i) = e->dot(*h);
79  }
80  }
81 
82  return H;
83 
84  }
85 
86 
87  template<class Real>
88  Teuchos::SerialDenseMatrix<int, Real> computeDotMatrix(const Vector<Real> &x) {
89 
90  int dim = x.dimension();
91  Teuchos::SerialDenseMatrix<int, Real> M(dim, dim);
92 
93  Teuchos::RCP<Vector<Real> > ei = x.clone();
94  Teuchos::RCP<Vector<Real> > ej = x.clone();
95 
96  for (int i=0; i<dim; i++) {
97  ei = x.basis(i);
98  for (int j=0; j<dim; j++) {
99  ej = x.basis(j);
100  M(j,i) = ej->dot(*ei);
101  }
102  }
103 
104  return M;
105 
106  }
107 
108  template<class Real>
109  std::vector<std::vector<Real> > computeEigenvalues(const Teuchos::SerialDenseMatrix<int, Real> & mat) {
110 
111  Teuchos::LAPACK<int, Real> lapack;
112 
113  Teuchos::SerialDenseMatrix<int, Real> mymat(Teuchos::Copy, mat);
114 
115  char jobvl = 'N';
116  char jobvr = 'N';
117 
118  int n = mat.numRows();
119 
120  std::vector<Real> real(n, 0);
121  std::vector<Real> imag(n, 0);
122  std::vector<std::vector<Real> > eigenvals;
123 
124  Real* vl = 0;
125  Real* vr = 0;
126 
127  int ldvl = 1;
128  int ldvr = 1;
129 
130  int lwork = 4*n;
131 
132  std::vector<Real> work(lwork, 0);
133 
134  int info = 0;
135 
136  lapack.GEEV(jobvl, jobvr, n, &mymat(0,0), n, &real[0], &imag[0], vl, ldvl, vr, ldvr, &work[0], lwork, &info);
137 
138  eigenvals.push_back(real);
139  eigenvals.push_back(imag);
140 
141  return eigenvals;
142 
143  }
144 
145 
146  template<class Real>
147  std::vector<std::vector<Real> > computeGenEigenvalues(const Teuchos::SerialDenseMatrix<int, Real> & A,
148  const Teuchos::SerialDenseMatrix<int, Real> & B) {
149 
150  Teuchos::LAPACK<int, Real> lapack;
151 
152  Teuchos::SerialDenseMatrix<int, Real> myA(Teuchos::Copy, A);
153  Teuchos::SerialDenseMatrix<int, Real> myB(Teuchos::Copy, B);
154 
155  char jobvl = 'N';
156  char jobvr = 'N';
157 
158  int n = A.numRows();
159 
160  std::vector<Real> real(n, 0);
161  std::vector<Real> imag(n, 0);
162  std::vector<Real> beta(n, 0);
163  std::vector<std::vector<Real> > eigenvals;
164 
165  Real* vl = 0;
166  Real* vr = 0;
167 
168  int ldvl = 1;
169  int ldvr = 1;
170 
171  int lwork = 10*n;
172 
173  std::vector<Real> work(lwork, 0);
174 
175  int info = 0;
176 
177  lapack.GGEV(jobvl, jobvr, n, &myA(0,0), n, &myB(0,0), n, &real[0], &imag[0], &beta[0],
178  vl, ldvl, vr, ldvr, &work[0], lwork, &info);
179 
180  for (int i=0; i<n; i++) {
181  real[i] /= beta[i];
182  imag[i] /= beta[i];
183  }
184 
185  eigenvals.push_back(real);
186  eigenvals.push_back(imag);
187 
188  return eigenvals;
189 
190  }
191 
192 
193  template<class Real>
194  Teuchos::SerialDenseMatrix<int, Real> computeInverse(const Teuchos::SerialDenseMatrix<int, Real> & mat) {
195 
196  Teuchos::LAPACK<int, Real> lapack;
197 
198  Teuchos::SerialDenseMatrix<int, Real> mymat(Teuchos::Copy, mat);
199 
200  int n = mat.numRows();
201 
202  std::vector<int> ipiv(n, 0);
203 
204  int lwork = 5*n;
205 
206  std::vector<Real> work(lwork, 0);
207 
208  int info = 0;
209 
210  lapack.GETRF(n, n, &mymat(0,0), n, &ipiv[0], &info);
211  lapack.GETRI(n, &mymat(0,0), n, &ipiv[0], &work[0], lwork, &info);
212 
213  return mymat;
214 
215  }
216 
217 
218 
219  template<class Real>
220  class ProjectedObjective : public Objective<Real> {
221  private:
222  Teuchos::RCP<Objective<Real> > obj_;
223  Teuchos::RCP<BoundConstraint<Real> > con_;
224  Teuchos::RCP<Secant<Real> > secant_;
227  Real eps_;
228 
229  public:
231  Teuchos::RCP<Secant<Real> > &secant,
232  bool useSecantPrecond = false, bool useSecantHessVec = false, Real eps = 0.0 ) {
233  obj_ = Teuchos::rcp(&obj, false);
234  con_ = Teuchos::rcp(&con, false);
235  secant_ = secant; //Teuchos::rcp(&secant, false);
236  useSecantPrecond_ = useSecantPrecond;
237  useSecantHessVec_ = useSecantHessVec;
238  eps_ = eps;
239  }
240 
241  void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
242  this->obj_->update(x,flag,iter);
243  this->con_->update(x,flag,iter);
244  }
245 
246  Real value( const Vector<Real> &x, Real &tol ) {
247  return this->obj_->value(x,tol);
248  }
249 
250  void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
251  this->obj_->gradient(g,x,tol);
252  }
253 
254  Real dirDeriv( const Vector<Real> &x, const Vector<Real> &d, Real &tol ) {
255  return this->obj_->dirDeriv(x,d,tol);
256  }
257 
258  void hessVec( Vector<Real> &Hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
259  if ( this->useSecantHessVec_ ) {
260  this->secant_->applyB( Hv, v, x );
261  }
262  else {
263  this->obj_->hessVec( Hv, v, x, tol );
264  }
265  }
266 
267  void invHessVec( Vector<Real> &Hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
268  if ( this->useSecantHessVec_ ) {
269  this->secant_->applyH(Hv,v,x);
270  }
271  else {
272  this->obj_->invHessVec(Hv,v,x,tol);
273  }
274  }
275 
276  void precond( Vector<Real> &Mv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
277  if ( this->useSecantPrecond_ ) {
278  this->secant_->applyH( Mv, v, x );
279  }
280  else {
281  this->obj_->precond( Mv, v, x, tol );
282  }
283  }
284 
296  void reducedHessVec( Vector<Real> &Hv, const Vector<Real> &v, const Vector<Real> &p,
297  const Vector<Real> &d, const Vector<Real> &x, Real &tol ) {
298  if ( this->con_->isActivated() ) {
299  Teuchos::RCP<Vector<Real> > vnew = x.clone();
300  // Set vnew to v
301  vnew->set(v);
302  // Remove elements of vnew corresponding to binding set
303  this->con_->pruneActive(*vnew,d,p,this->eps_);
304  // Apply full Hessian to reduced vector
305  this->hessVec(Hv,*vnew,x,tol);
306  // Remove elements of Hv corresponding to binding set
307  this->con_->pruneActive(Hv,d,p,this->eps_);
308  // Set vnew to v
309  vnew->set(v);
310  // Remove Elements of vnew corresponding to complement of binding set
311  this->con_->pruneInactive(*vnew,d,p,this->eps_);
312  // Fill complement of binding set elements in Hp with v
313  Hv.plus(*vnew);
314  }
315  else {
316  this->hessVec(Hv,v,x,tol);
317  }
318  }
319 
330  void reducedHessVec( Vector<Real> &Hv, const Vector<Real> &v, const Vector<Real> &p,
331  const Vector<Real> &x, Real &tol ) {
332  if ( this->con_->isActivated() ) {
333  Teuchos::RCP<Vector<Real> > vnew = x.clone();
334  // Set vnew to v
335  vnew->set(v);
336  // Remove elements of vnew corresponding to binding set
337  this->con_->pruneActive(*vnew,p,this->eps_);
338  // Apply full Hessian to reduced vector
339  this->hessVec(Hv,*vnew,x,tol);
340  // Remove elements of Hv corresponding to binding set
341  this->con_->pruneActive(Hv,p,this->eps_);
342  // Set vnew to v
343  vnew->set(v);
344  // Remove Elements of vnew corresponding to complement of binding set
345  this->con_->pruneInactive(*vnew,p,this->eps_);
346  // Fill complement of binding set elements in Hp with v
347  Hv.plus(*vnew);
348  }
349  else {
350  this->hessVec(Hv,v,x,tol);
351  }
352  }
353 
365  void reducedInvHessVec( Vector<Real> &Hv, const Vector<Real> &v, const Vector<Real> &p,
366  const Vector<Real> &d, const Vector<Real> &x, Real &tol ) {
367  if ( this->con_->isActivated() ) {
368  Teuchos::RCP<Vector<Real> > vnew = x.clone();
369  // Set vnew to v
370  vnew->set(v);
371  // Remove elements of vnew corresponding to binding set
372  this->con_->pruneActive(*vnew,d,p,this->eps_);
373  // Apply full Hessian to reduced vector
374  this->invHessVec(Hv,*vnew,x,tol);
375  // Remove elements of Hv corresponding to binding set
376  this->con_->pruneActive(Hv,d,p,this->eps_);
377  // Set vnew to v
378  vnew->set(v);
379  // Remove Elements of vnew corresponding to complement of binding set
380  this->con_->pruneInactive(*vnew,d,p,this->eps_);
381  // Fill complement of binding set elements in Hv with v
382  Hv.plus(*vnew);
383  }
384  else {
385  this->invHessVec(Hv,v,x,tol);
386  }
387  }
388 
399  void reducedInvHessVec( Vector<Real> &Hv, const Vector<Real> &v, const Vector<Real> &p,
400  const Vector<Real> &x, Real &tol ) {
401  if ( this->con_->isActivated() ) {
402  Teuchos::RCP<Vector<Real> > vnew = x.clone();
403  // Set vnew to v
404  vnew->set(v);
405  // Remove elements of vnew corresponding to binding set
406  this->con_->pruneActive(*vnew,p,this->eps_);
407  // Apply full Hessian to reduced vector
408  this->invHessVec(Hv,*vnew,x,tol);
409  // Remove elements of Hv corresponding to binding set
410  this->con_->pruneActive(Hv,p,this->eps_);
411  // Set vnew to v
412  vnew->set(v);
413  // Remove Elements of vnew corresponding to complement of binding set
414  this->con_->pruneInactive(*vnew,p,this->eps_);
415  // Fill complement of binding set elements in Hv with v
416  Hv.plus(*vnew);
417  }
418  else {
419  this->invHessVec(Hv,v,x,tol);
420  }
421  }
422 
434  void reducedPrecond( Vector<Real> &Mv, const Vector<Real> &v, const Vector<Real> &p,
435  const Vector<Real> &d, const Vector<Real> &x, Real &tol ) {
436  if ( this->con_->isActivated() ) {
437  Teuchos::RCP<Vector<Real> > vnew = x.clone();
438  // Set vnew to v
439  vnew->set(v);
440  // Remove elements of vnew corresponding to binding set
441  this->con_->pruneActive(*vnew,d,p,this->eps_);
442  // Apply full Hessian to reduced vector
443  this->precond(Mv,*vnew,x,tol);
444  // Remove elements of Mv corresponding to binding set
445  this->con_->pruneActive(Mv,d,p,this->eps_);
446  // Set vnew to v
447  vnew->set(v);
448  // Remove Elements of vnew corresponding to complement of binding set
449  this->con_->pruneInactive(*vnew,d,p,this->eps_);
450  // Fill complement of binding set elements in Mv with v
451  Mv.plus(*vnew);
452  }
453  else {
454  this->precond(Mv,v,x,tol);
455  }
456  }
457 
468  void reducedPrecond( Vector<Real> &Mv, const Vector<Real> &v, const Vector<Real> &p,
469  const Vector<Real> &x, Real &tol ) {
470  if ( this->con_->isActivated() ) {
471  Teuchos::RCP<Vector<Real> > vnew = x.clone();
472  // Set vnew to v
473  vnew->set(v);
474  // Remove elements of vnew corresponding to binding set
475  this->con_->pruneActive(*vnew,p,this->eps_);
476  // Apply full Hessian to reduced vector
477  this->precond(Mv,*vnew,x,tol);
478  // Remove elements of Mv corresponding to binding set
479  this->con_->pruneActive(Mv,p,this->eps_);
480  // Set vnew to v
481  vnew->set(v);
482  // Remove Elements of vnew corresponding to complement of binding set
483  this->con_->pruneInactive(*vnew,p,this->eps_);
484  // Fill complement of binding set elements in Mv with v
485  Mv.plus(*vnew);
486  }
487  else {
488  this->precond(Mv,v,x,tol);
489  }
490  }
491 
492  void project( Vector<Real> &x ) {
493  this->con_->project(x);
494  }
495 
496  void pruneActive( Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x ) {
497  this->con_->pruneActive(v,g,x,this->eps_);
498  }
499 
500  void pruneActive( Vector<Real> &v, const Vector<Real> &x ) {
501  this->con_->pruneActive(v,x,this->eps_);
502  }
503 
504  void pruneInactive( Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x ) {
505  this->con_->pruneInactive(v,g,x,this->eps_);
506  }
507 
508  void pruneInactive( Vector<Real> &v, const Vector<Real> &x ) {
509  this->con_->pruneInactive(v,x,this->eps_);
510  }
511 
512  bool isFeasible( const Vector<Real> &v ) {
513  return this->con_->isFeasible(v);
514  }
515 
516  bool isConActivated(void) {
517  return this->con_->isActivated();
518  }
519 
521  this->con_->computeProjectedStep(v,x);
522  }
523  };
524 
525 } // namespace ROL
526 
527 #endif
Provides the interface to evaluate objective functions.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
virtual int dimension() const
Return dimension of the vector space.
Definition: ROL_Vector.hpp:181
virtual void plus(const Vector &x)=0
Compute , where .
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...
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Real dirDeriv(const Vector< Real > &x, const Vector< Real > &d, Real &tol)
Compute directional derivative.
void pruneInactive(Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x)
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
Teuchos::RCP< BoundConstraint< Real > > con_
bool isFeasible(const Vector< Real > &v)
void pruneActive(Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x)
void reducedPrecond(Vector< Real > &Mv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &x, Real &tol)
Apply the reduced preconditioner to a vector, v. The reduced preconditioner first removes elements of...
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
void hessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Teuchos::SerialDenseMatrix< int, Real > computeDenseHessian(Objective< Real > &obj, const Vector< Real > &x)
std::vector< std::vector< Real > > computeGenEigenvalues(const Teuchos::SerialDenseMatrix< int, Real > &A, const Teuchos::SerialDenseMatrix< int, Real > &B)
void reducedInvHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &x, Real &tol)
Apply the reduced inverse Hessian to a vector, v. The reduced inverse Hessian first removes elements ...
ProjectedObjective(Objective< Real > &obj, BoundConstraint< Real > &con, Teuchos::RCP< Secant< Real > > &secant, bool useSecantPrecond=false, bool useSecantHessVec=false, Real eps=0.0)
Provides interface for and implements limited-memory secant operators.
Definition: ROL_Secant.hpp:67
Teuchos::RCP< Objective< Real > > obj_
std::vector< std::vector< Real > > computeEigenvalues(const Teuchos::SerialDenseMatrix< int, Real > &mat)
void reducedInvHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &d, const Vector< Real > &x, Real &tol)
Apply the reduced inverse Hessian to a vector, v. The reduced inverse Hessian first removes elements ...
void pruneInactive(Vector< Real > &v, const Vector< Real > &x)
void reducedHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &x, Real &tol)
Apply the reduced Hessian to a vector, v. The reduced Hessian first removes elements of v correspondi...
Provides the interface to apply upper and lower bound constraints.
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
Teuchos::SerialDenseMatrix< int, Real > computeDotMatrix(const Vector< Real > &x)
virtual Teuchos::RCP< Vector > basis(const int i) const
Return i-th basis vector.
Definition: ROL_Vector.hpp:170
void computeProjectedStep(Vector< Real > &v, const Vector< Real > &x)
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...
void project(Vector< Real > &x)
void invHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
Teuchos::SerialDenseMatrix< int, Real > computeInverse(const Teuchos::SerialDenseMatrix< int, Real > &mat)
void precond(Vector< Real > &Mv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply preconditioner to vector.
void pruneActive(Vector< Real > &v, const Vector< Real > &x)
Teuchos::RCP< Secant< Real > > secant_
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:115