ROL
ROL_GMRES.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_GMRES_H
11#define ROL_GMRES_H
12
17#include "ROL_Krylov.hpp"
18#include "ROL_Types.hpp"
19#include "ROL_LAPACK.hpp"
20#include "ROL_LinearAlgebra.hpp"
21
22
23namespace ROL {
24
25template<class Real>
27private:
28 std::vector<Ptr<Vector<Real>>> data_;
29 unsigned size_;
30
31public:
32 VectorArray() : size_(0u) {}
33
34 void allocate(const Vector<Real> &x, unsigned numVec) {
35 data_.clear();
36 for (unsigned i = 0u; i < numVec; ++i)
37 data_.push_back(x.clone());
38 size_ = numVec;
39 }
40
41 const Ptr<Vector<Real>> get(const Vector<Real> &x, unsigned ind) {
42 if (ind < size_) {
43 return data_[ind];
44 }
45 else if (ind == size_) {
46 data_.push_back(x.clone());
47 size_++;
48 return data_[ind];
49 }
50 else {
51 throw Exception::NotImplemented(">>> GMRES : Index out of range!");
52 return nullPtr;
53 }
54 }
55
56 const Ptr<Vector<Real>> get(unsigned ind) {
57 if (size_ > 0u) {
58 return get(*data_[0],ind);
59 }
60 else {
61 throw Exception::NotImplemented(">>> GMRES : No vectors allocated!");
62 return nullPtr;
63 }
64 }
65};
66
67template<class Real>
68class GMRES : public Krylov<Real> {
69
70 typedef LA::Matrix<Real> SDMatrix;
71 typedef LA::Vector<Real> SDVector;
72
73private:
74
75 Ptr<Vector<Real>> r_;
76 Ptr<Vector<Real>> z_;
77 Ptr<Vector<Real>> w_;
78
79 Ptr<SDMatrix> H_; // quasi-Hessenberg matrix
80 Ptr<SDVector> cs_; // Givens Rotations cosine components
81 Ptr<SDVector> sn_; // Givens Rotations sine components
82 Ptr<SDVector> s_;
83 Ptr<SDVector> y_;
84 Ptr<SDVector> cnorm_;
85
86 Ptr<std::vector<Real>> res_;
87 Ptr<VectorArray<Real>> V_, Z_;
88
91 bool useInitialGuess_; // If false, inital x will be ignored and zero vec used
93 Ptr<std::ostream> outStream_;
94
95 LAPACK<int,Real> lapack_;
96
97public:
98
99 GMRES( ParameterList &parlist ) : Krylov<Real>(parlist), isInitialized_(false), printIters_(false) {
100
101 using std::vector;
102
103 Real zero(0);
104
105 ROL::ParameterList &gList = parlist.sublist("General");
106 ROL::ParameterList &kList = gList.sublist("Krylov");
107
108 useInexact_ = gList.get("Inexact Hessian-Times-A-Vector",false);
109 useInitialGuess_ = kList.get("Use Initial Guess",false);
111
112 H_ = makePtr<SDMatrix>( maxit+1, maxit );
113 cs_ = makePtr<SDVector>( maxit );
114 sn_ = makePtr<SDVector>( maxit );
115 s_ = makePtr<SDVector>( maxit+1 );
116 y_ = makePtr<SDVector>( maxit+1 );
117 cnorm_ = makePtr<SDVector>( maxit );
118 res_ = makePtr<std::vector<Real>>(maxit+1,zero);
119
120 V_ = makePtr<VectorArray<Real>>();
121 Z_ = makePtr<VectorArray<Real>>();
122 }
123
125 LinearOperator<Real> &M, int &iter, int &flag ) {
126
130
131 flag = 0;
132
133 Real zero(0), one(1);
134
135 if ( !isInitialized_ ) {
136 r_ = b.clone();
137 w_ = b.clone();
138 z_ = x.clone();
139
140 V_->allocate(b,std::min(20u,static_cast<unsigned>(maxit))+1u);
141 Z_->allocate(x,std::min(20u,static_cast<unsigned>(maxit)));
142
143 isInitialized_ = true;
144 }
145
146 Real itol = std::sqrt(ROL_EPSILON<Real>());
147
148 // Compute initial residual
149 if(useInitialGuess_) {
150
151 A.apply(*r_,x,itol);
152 r_->scale(-one);
153 r_->plus(b); // r = b-Ax
154
155 }
156 else {
157 x.zero();
158 r_->set(b);
159 }
160
161 Real temp = 0;
162
163 (*res_)[0] = r_->norm();
164
165 if (printIters_)
166 *outStream_ << "GMRES Iteration " << 0 << ", Residual = " << (*res_)[0] << "\n";
167
168 // This should be a tolerance check
169 Real rtol = std::min(absTol,relTol*(*res_)[0]);
170 if ((*res_)[0] <= rtol) {
171 iter = 0;
172 flag = 0;
173 return (*res_)[0];
174 }
175
176 V_->get(0)->set(*r_);
177 V_->get(0)->scale(one/(*res_)[0]);
178
179 (*s_)(0) = (*res_)[0];
180
181 for( iter=0; iter<maxit; ++iter ) {
182
183 if( useInexact_ ) itol = rtol/(maxit*(*res_)[iter]);
184
185 // Apply right preconditioner
186 M.applyInverse(*(Z_->get(iter)),*(V_->get(iter)),itol);
187
188 // Apply operator
189 A.apply(*w_,*(Z_->get(iter)),itol);
190
191 // Evaluate coefficients and orthogonalize using Gram-Schmidt
192 for( int k=0; k<=iter; ++k ) {
193 (*H_)(k,iter) = w_->dot(*(V_->get(k)));
194 w_->axpy( -(*H_)(k,iter), *(V_->get(k)) );
195 }
196
197 (*H_)(iter+1,iter) = w_->norm();
198
199 (V_->get(iter+1))->set(*w_);
200 (V_->get(iter+1))->scale(one/((*H_)(iter+1,iter)));
201
202 // Apply Givens rotations
203 for( int k=0; k<=iter-1; ++k ) {
204 temp = (*cs_)(k)*(*H_)(k,iter) + (*sn_)(k)*(*H_)(k+1,iter);
205 (*H_)(k+1,iter) = -(*sn_)(k)*(*H_)(k,iter) + (*cs_)(k)*(*H_)(k+1,iter);
206 (*H_)(k,iter) = temp;
207 }
208
209 // Form i-th rotation matrix
210 if( (*H_)(iter+1,iter) == zero ) {
211 (*cs_)(iter) = one;
212 (*sn_)(iter) = zero;
213 }
214 else if ( std::abs((*H_)(iter+1,iter)) > std::abs((*H_)(iter,iter)) ) {
215 temp = (*H_)(iter,iter) / (*H_)(iter+1,iter);
216 (*sn_)(iter) = one / std::sqrt( one + temp*temp );
217 (*cs_)(iter) = temp*(*sn_)(iter);
218 }
219 else {
220 temp = (*H_)(iter+1,iter) / (*H_)(iter,iter);
221 (*cs_)(iter) = one / std::sqrt( one + temp*temp );
222 (*sn_)(iter) = temp*(*cs_)(iter);
223 }
224
225 // Approximate residual norm
226 temp = (*cs_)(iter)*(*s_)(iter);
227 (*s_)(iter+1) = -(*sn_)(iter)*(*s_)(iter);
228 (*s_)(iter) = temp;
229 (*H_)(iter,iter) = (*cs_)(iter)*(*H_)(iter,iter) + (*sn_)(iter)*(*H_)(iter+1,iter);
230 (*H_)(iter+1,iter) = zero;
231 (*res_)[iter+1] = std::abs((*s_)(iter+1));
232
233 if (printIters_) {
234 *outStream_ << "GMRES Iteration " << iter+1 << ", Residual = " << (*res_)[iter+1] << "\n";
235 }
236
237 // Update solution approximation.
238 const char uplo = 'U';
239 const char trans = 'N';
240 const char diag = 'N';
241 const char normin = 'N';
242 Real scaling = zero;
243 int info = 0;
244 *y_ = *s_;
245 lapack_.LATRS(uplo, trans, diag, normin, iter+1, H_->values(), maxit+1, y_->values(), &scaling, cnorm_->values(), &info);
246
247 z_->zero();
248
249 for( int k=0; k<=iter;++k ) {
250 z_->axpy((*y_)(k),*(Z_->get(k)));
251 }
252
253 if( (*res_)[iter+1] <= rtol ) {
254 // Update solution vector
255 x.plus(*z_);
256 break;
257 }
258
259 } // loop over iter
260
261 if(iter == maxit) {
262 flag = 1;
263 x.plus(*z_);
264 return (*res_)[iter];
265 }
266
267 return (*res_)[iter+1];
268 }
269
270 void enableOutput(std::ostream & outStream) {
271 printIters_ = true;
272 outStream_ = ROL::makePtrFromRef(outStream);;
273 }
274
275 void disableOutput() {printIters_ = false;}
276
277}; // class GMRES
278
279} // namespace ROL
280
281#endif // ROL_GMRES_H
282
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Contains definitions of custom data types in ROL.
Preconditioned GMRES solver.
Definition ROL_GMRES.hpp:68
GMRES(ParameterList &parlist)
Definition ROL_GMRES.hpp:99
LA::Vector< Real > SDVector
Definition ROL_GMRES.hpp:71
Ptr< VectorArray< Real > > Z_
Definition ROL_GMRES.hpp:87
LAPACK< int, Real > lapack_
Definition ROL_GMRES.hpp:95
Ptr< SDMatrix > H_
Definition ROL_GMRES.hpp:79
Ptr< SDVector > sn_
Definition ROL_GMRES.hpp:81
Ptr< SDVector > cnorm_
Definition ROL_GMRES.hpp:84
Ptr< Vector< Real > > r_
Definition ROL_GMRES.hpp:75
void enableOutput(std::ostream &outStream)
Ptr< Vector< Real > > w_
Definition ROL_GMRES.hpp:77
Ptr< SDVector > s_
Definition ROL_GMRES.hpp:82
bool useInitialGuess_
Definition ROL_GMRES.hpp:91
bool useInexact_
Definition ROL_GMRES.hpp:90
Ptr< std::ostream > outStream_
Definition ROL_GMRES.hpp:93
void disableOutput()
LA::Matrix< Real > SDMatrix
Definition ROL_GMRES.hpp:70
Ptr< Vector< Real > > z_
Definition ROL_GMRES.hpp:76
bool isInitialized_
Definition ROL_GMRES.hpp:89
Ptr< std::vector< Real > > res_
Definition ROL_GMRES.hpp:86
Ptr< VectorArray< Real > > V_
Definition ROL_GMRES.hpp:87
Ptr< SDVector > y_
Definition ROL_GMRES.hpp:83
Real run(Vector< Real > &x, LinearOperator< Real > &A, const Vector< Real > &b, LinearOperator< Real > &M, int &iter, int &flag)
Ptr< SDVector > cs_
Definition ROL_GMRES.hpp:80
bool printIters_
Definition ROL_GMRES.hpp:92
Provides definitions for Krylov solvers.
unsigned getMaximumIteration(void) const
Real getAbsoluteTolerance(void) const
Real getRelativeTolerance(void) const
Provides the interface to apply a linear operator.
virtual void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const =0
Apply linear operator.
virtual void applyInverse(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply inverse of linear operator.
std::vector< Ptr< Vector< Real > > > data_
Definition ROL_GMRES.hpp:28
const Ptr< Vector< Real > > get(const Vector< Real > &x, unsigned ind)
Definition ROL_GMRES.hpp:41
void allocate(const Vector< Real > &x, unsigned numVec)
Definition ROL_GMRES.hpp:34
const Ptr< Vector< Real > > get(unsigned ind)
Definition ROL_GMRES.hpp:56
Defines the linear algebra or vector space interface.
virtual void plus(const Vector &x)=0
Compute , where .
virtual void zero()
Set to zero vector.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.