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 47 /** \class ROL::GMRES 48 \brief Preconditioned GMRES solver. 49 */ 50 51 #include "ROL_Krylov.hpp" 52 #include "ROL_Types.hpp" 53 #include "ROL_LAPACK.hpp" 54 #include "ROL_LinearAlgebra.hpp" 55 56 57 namespace ROL { 58 59 template<class Real> 60 class GMRES : public Krylov<Real> { 61 62 typedef LA::Matrix<Real> SDMatrix; 63 typedef LA::Vector<Real> SDVector; 64 65 private: 66 67 ROL::Ptr<Vector<Real> > r_; 68 ROL::Ptr<Vector<Real> > z_; 69 ROL::Ptr<Vector<Real> > w_; 70 71 ROL::Ptr<SDMatrix> H_; // quasi-Hessenberg matrix 72 ROL::Ptr<SDVector> cs_; // Givens Rotations cosine components 73 ROL::Ptr<SDVector> sn_; // Givens Rotations sine components 74 ROL::Ptr<SDVector> s_; 75 ROL::Ptr<SDVector> y_; 76 ROL::Ptr<SDVector> cnorm_; 77 78 ROL::Ptr<std::vector<Real> > res_; 79 80 bool isInitialized_; 81 bool useInexact_; 82 bool useInitialGuess_; // If false, inital x will be ignored and zero vec used 83 bool printIters_; 84 ROL::Ptr<std::ostream> outStream_; 85 86 ROL::LAPACK<int,Real> lapack_; 87 88 public: 89 GMRES(ROL::ParameterList & parlist)90 GMRES( ROL::ParameterList &parlist ) : Krylov<Real>(parlist), isInitialized_(false), printIters_(false) { 91 92 using std::vector; 93 94 Real zero(0); 95 96 ROL::ParameterList &gList = parlist.sublist("General"); 97 ROL::ParameterList &kList = gList.sublist("Krylov"); 98 99 useInexact_ = gList.get("Inexact Hessian-Times-A-Vector",false); 100 useInitialGuess_ = kList.get("Use Initial Guess",false); 101 int maxit = Krylov<Real>::getMaximumIteration(); 102 103 H_ = ROL::makePtr<SDMatrix>( maxit+1, maxit ); 104 cs_ = ROL::makePtr<SDVector>( maxit ); 105 sn_ = ROL::makePtr<SDVector>( maxit ); 106 s_ = ROL::makePtr<SDVector>( maxit+1 ); 107 y_ = ROL::makePtr<SDVector>( maxit+1 ); 108 cnorm_ = ROL::makePtr<SDVector>( maxit ); 109 res_ = ROL::makePtr<std::vector<Real>>(maxit+1,zero); 110 111 } 112 run(Vector<Real> & x,LinearOperator<Real> & A,const Vector<Real> & b,LinearOperator<Real> & M,int & iter,int & flag)113 Real run( Vector<Real> &x, LinearOperator<Real> &A, const Vector<Real> &b, 114 LinearOperator<Real> &M, int &iter, int &flag ) { 115 116 Real absTol = Krylov<Real>::getAbsoluteTolerance(); 117 Real relTol = Krylov<Real>::getRelativeTolerance(); 118 int maxit = Krylov<Real>::getMaximumIteration(); 119 120 flag = 0; 121 122 Real zero(0), one(1); 123 124 if ( !isInitialized_ ) { 125 r_ = b.clone(); 126 w_ = b.clone(); 127 z_ = x.clone(); 128 129 isInitialized_ = true; 130 } 131 132 Real itol = std::sqrt(ROL_EPSILON<Real>()); 133 134 // Compute initial residual 135 if(useInitialGuess_) { 136 137 A.apply(*r_,x,itol); 138 r_->scale(-one); 139 r_->plus(b); // r = b-Ax 140 141 } 142 else { 143 x.zero(); 144 r_->set(b); 145 } 146 147 Real temp = 0; 148 149 std::vector<ROL::Ptr<Vector<Real > > > V; 150 std::vector<ROL::Ptr<Vector<Real > > > Z; 151 152 (*res_)[0] = r_->norm(); 153 154 if (printIters_) { 155 *outStream_ << "GMRES Iteration " << 0 << ", Residual = " << (*res_)[0] << "\n"; 156 } 157 158 // This should be a tolerance check 159 Real rtol = std::min(absTol,relTol*(*res_)[0]); 160 if ((*res_)[0] <= rtol) { 161 iter = 0; 162 flag = 0; 163 return (*res_)[0]; 164 } 165 166 V.push_back(b.clone()); 167 (V[0])->set(*r_); 168 (V[0])->scale(one/(*res_)[0]); 169 170 (*s_)(0) = (*res_)[0]; 171 172 for( iter=0; iter<maxit; ++iter ) { 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 if (printIters_) { 230 *outStream_ << "GMRES Iteration " << iter+1 << ", Residual = " << (*res_)[iter+1] << "\n"; 231 } 232 233 // Update solution approximation. 234 const char uplo = 'U'; 235 const char trans = 'N'; 236 const char diag = 'N'; 237 const char normin = 'N'; 238 Real scaling = zero; 239 int info = 0; 240 *y_ = *s_; 241 lapack_.LATRS(uplo, trans, diag, normin, iter+1, H_->values(), maxit+1, y_->values(), &scaling, cnorm_->values(), &info); 242 243 z_->zero(); 244 245 for( int k=0; k<=iter;++k ) { 246 z_->axpy((*y_)(k),*(Z[k])); 247 } 248 249 if( (*res_)[iter+1] <= rtol ) { 250 // Update solution vector 251 x.plus(*z_); 252 break; 253 } 254 255 } // loop over iter 256 257 if(iter == maxit) { 258 flag = 1; 259 x.plus(*z_); 260 return (*res_)[iter]; 261 } 262 263 return (*res_)[iter+1]; 264 } 265 enableOutput(std::ostream & outStream)266 void enableOutput(std::ostream & outStream) { 267 printIters_ = true; 268 outStream_ = ROL::makePtrFromRef(outStream);; 269 } 270 disableOutput()271 void disableOutput() {printIters_ = false;} 272 273 }; // class GMRES 274 275 } // namespace ROL 276 277 #endif // ROL_GMRES_H 278 279