ROL
ROL_GMRES.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_GMRES_H
45 #define ROL_GMRES_H
46 
51 #include "ROL_Krylov.hpp"
52 #include "ROL_LinearOperator.hpp"
53 #include "ROL_Vector.hpp"
54 #include "ROL_Types.hpp"
55 #include "ROL_LAPACK.hpp"
56 #include "ROL_LinearAlgebra.hpp"
57 
58 
59 namespace ROL {
60 
61 template<class Real>
62 class GMRES : public Krylov<Real> {
63 
64  typedef LA::Matrix<Real> SDMatrix;
65  typedef LA::Vector<Real> SDVector;
66 
67 private:
68 
69  ROL::Ptr<Vector<Real> > r_;
70  ROL::Ptr<Vector<Real> > z_;
71  ROL::Ptr<Vector<Real> > w_;
72 
73  ROL::Ptr<SDMatrix> H_; // quasi-Hessenberg matrix
74  ROL::Ptr<SDVector> cs_; // Givens Rotations cosine components
75  ROL::Ptr<SDVector> sn_; // Givens Rotations sine components
76  ROL::Ptr<SDVector> s_;
77  ROL::Ptr<SDVector> y_;
78  ROL::Ptr<SDVector> cnorm_;
79 
80  ROL::Ptr<std::vector<Real> > res_;
81 
84  bool useInitialGuess_; // If false, inital x will be ignored and zero vec used
85 
86  ROL::LAPACK<int,Real> lapack_;
87 
88 public:
89 
90  GMRES( ROL::ParameterList &parlist ) : Krylov<Real>(parlist), isInitialized_(false) {
91 
92 
93 
94  using std::vector;
95 
96  Real zero(0);
97 
98  ROL::ParameterList &gList = parlist.sublist("General");
99  ROL::ParameterList &kList = gList.sublist("Krylov");
100 
101  useInexact_ = gList.get("Inexact Hessian-Times-A-Vector",false);
102  useInitialGuess_ = kList.get("Use Initial Guess",false);
103  int maxit = Krylov<Real>::getMaximumIteration();
104 
105  H_ = ROL::makePtr<SDMatrix>( maxit+1, maxit );
106  cs_ = ROL::makePtr<SDVector>( maxit );
107  sn_ = ROL::makePtr<SDVector>( maxit );
108  s_ = ROL::makePtr<SDVector>( maxit+1 );
109  y_ = ROL::makePtr<SDVector>( maxit+1 );
110  cnorm_ = ROL::makePtr<SDVector>( maxit );
111  res_ = ROL::makePtr<std::vector<Real>>(maxit+1,zero);
112 
113  }
114 
116  LinearOperator<Real> &M, int &iter, int &flag ) {
117 
118  Real absTol = Krylov<Real>::getAbsoluteTolerance();
119  Real relTol = Krylov<Real>::getRelativeTolerance();
120  int maxit = Krylov<Real>::getMaximumIteration();
121 
122  flag = 0;
123 
124  Real zero(0), one(1);
125 
126  if ( !isInitialized_ ) {
127  r_ = b.clone();
128  w_ = b.clone();
129  z_ = x.clone();
130 
131  isInitialized_ = true;
132  }
133 
134  Real itol = std::sqrt(ROL_EPSILON<Real>());
135 
136  // Compute initial residual
137  if(useInitialGuess_) {
138 
139  A.apply(*r_,x,itol);
140  r_->scale(-one);
141  r_->plus(b); // r = b-Ax
142 
143  }
144  else {
145  x.zero();
146  r_->set(b);
147  }
148 
149  Real temp = 0;
150 
151  std::vector<ROL::Ptr<Vector<Real > > > V;
152  std::vector<ROL::Ptr<Vector<Real > > > Z;
153 
154  (*res_)[0] = r_->norm();
155 
156  // This should be a tolerance check
157  Real rtol = std::min(absTol,relTol*(*res_)[0]);
158  if ((*res_)[0] <= rtol) {
159  iter = 0;
160  flag = 0;
161  return (*res_)[0];
162  }
163 
164  V.push_back(b.clone());
165  (V[0])->set(*r_);
166  (V[0])->scale(one/(*res_)[0]);
167 
168  (*s_)(0) = (*res_)[0];
169 
170  for( iter=0; iter<maxit; ++iter ) {
171 
172  // std::cout << "iter = " << iter << " rnorm = " << (*res_)[iter] << " xnorm = " << x.norm() << " bnorm = " << b.norm() << std::endl;
173 
174  if( useInexact_ ) {
175  itol = rtol/(maxit*(*res_)[iter]);
176  }
177 
178  Z.push_back(x.clone());
179 
180  // Apply right preconditioner
181  M.applyInverse(*(Z[iter]),*(V[iter]),itol);
182 
183  // Apply operator
184  A.apply(*w_,*(Z[iter]),itol);
185 
186  // Evaluate coefficients and orthogonalize using Gram-Schmidt
187  for( int k=0; k<=iter; ++k ) {
188  (*H_)(k,iter) = w_->dot(*(V[k]));
189  w_->axpy( -(*H_)(k,iter), *(V[k]) );
190  }
191 
192  (*H_)(iter+1,iter) = w_->norm();
193 
194  V.push_back( b.clone() );
195  (V[iter+1])->set(*w_);
196  (V[iter+1])->scale(one/((*H_)(iter+1,iter)));
197 
198  // Apply Givens rotations
199  for( int k=0; k<=iter-1; ++k ) {
200  temp = (*cs_)(k)*(*H_)(k,iter) + (*sn_)(k)*(*H_)(k+1,iter);
201  (*H_)(k+1,iter) = -(*sn_)(k)*(*H_)(k,iter) + (*cs_)(k)*(*H_)(k+1,iter);
202  (*H_)(k,iter) = temp;
203  }
204 
205  // Form i-th rotation matrix
206  if( (*H_)(iter+1,iter) == zero ) {
207  (*cs_)(iter) = one;
208  (*sn_)(iter) = zero;
209  }
210  else if ( std::abs((*H_)(iter+1,iter)) > std::abs((*H_)(iter,iter)) ) {
211  temp = (*H_)(iter,iter) / (*H_)(iter+1,iter);
212  (*sn_)(iter) = one / std::sqrt( one + temp*temp );
213  (*cs_)(iter) = temp*(*sn_)(iter);
214  }
215  else {
216  temp = (*H_)(iter+1,iter) / (*H_)(iter,iter);
217  (*cs_)(iter) = one / std::sqrt( one + temp*temp );
218  (*sn_)(iter) = temp*(*cs_)(iter);
219  }
220 
221  // Approximate residual norm
222  temp = (*cs_)(iter)*(*s_)(iter);
223  (*s_)(iter+1) = -(*sn_)(iter)*(*s_)(iter);
224  (*s_)(iter) = temp;
225  (*H_)(iter,iter) = (*cs_)(iter)*(*H_)(iter,iter) + (*sn_)(iter)*(*H_)(iter+1,iter);
226  (*H_)(iter+1,iter) = zero;
227  (*res_)[iter+1] = std::abs((*s_)(iter+1));
228 
229  // Update solution approximation.
230  const char uplo = 'U';
231  const char trans = 'N';
232  const char diag = 'N';
233  const char normin = 'N';
234  Real scaling = zero;
235  int info = 0;
236  *y_ = *s_;
237  lapack_.LATRS(uplo, trans, diag, normin, iter+1, H_->values(), maxit+1, y_->values(), &scaling, cnorm_->values(), &info);
238 
239  z_->zero();
240 
241  for( int k=0; k<=iter;++k ) {
242  z_->axpy((*y_)(k),*(Z[k]));
243  }
244 
245  if( (*res_)[iter+1] <= rtol ) {
246  // Update solution vector
247  x.plus(*z_);
248  break;
249  }
250  } // loop over iter
251 
252  if(iter == maxit) {
253  flag = 1;
254  x.plus(*z_);
255  return (*res_)[iter];
256  }
257 
258  return (*res_)[iter+1];
259  }
260 
261 }; // class GMRES
262 
263 } // namespace ROL
264 
265 #endif // ROL_GMRES_H
266 
ROL::Ptr< SDVector > cs_
Definition: ROL_GMRES.hpp:74
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
Real getRelativeTolerance(void) const
Definition: ROL_Krylov.hpp:93
Real getAbsoluteTolerance(void) const
Definition: ROL_Krylov.hpp:90
bool isInitialized_
Definition: ROL_GMRES.hpp:82
virtual void plus(const Vector &x)=0
Compute , where .
GMRES(ROL::ParameterList &parlist)
Definition: ROL_GMRES.hpp:90
bool useInexact_
Definition: ROL_GMRES.hpp:83
Contains definitions of custom data types in ROL.
virtual void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const =0
Apply linear operator.
unsigned getMaximumIteration(void) const
Definition: ROL_Krylov.hpp:96
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:167
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
ROL::Ptr< Vector< Real > > r_
Definition: ROL_GMRES.hpp:69
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Vector< Real > V
LA::Matrix< Real > SDMatrix
Definition: ROL_GMRES.hpp:64
Real run(Vector< Real > &x, LinearOperator< Real > &A, const Vector< Real > &b, LinearOperator< Real > &M, int &iter, int &flag)
Definition: ROL_GMRES.hpp:115
ROL::Ptr< std::vector< Real > > res_
Definition: ROL_GMRES.hpp:80
Preconditioned GMRES solver.
Definition: ROL_GMRES.hpp:62
bool useInitialGuess_
Definition: ROL_GMRES.hpp:84
virtual void applyInverse(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply inverse of linear operator.
ROL::Ptr< SDVector > y_
Definition: ROL_GMRES.hpp:77
Provides definitions for Krylov solvers.
Definition: ROL_Krylov.hpp:58
Provides the interface to apply a linear operator.
ROL::Ptr< SDVector > cnorm_
Definition: ROL_GMRES.hpp:78
ROL::Ptr< SDVector > sn_
Definition: ROL_GMRES.hpp:75
ROL::Ptr< Vector< Real > > w_
Definition: ROL_GMRES.hpp:71
ROL::Ptr< SDMatrix > H_
Definition: ROL_GMRES.hpp:73
LA::Vector< Real > SDVector
Definition: ROL_GMRES.hpp:65
ROL::Ptr< Vector< Real > > z_
Definition: ROL_GMRES.hpp:70
ROL::Ptr< SDVector > s_
Definition: ROL_GMRES.hpp:76
ROL::LAPACK< int, Real > lapack_
Definition: ROL_GMRES.hpp:86