lanczos.h

00001 /*
00002   -------------------------------------------------------------------
00003   
00004   Copyright (C) 2006, 2007, 2008, Andrew W. Steiner
00005   
00006   This file is part of O2scl.
00007   
00008   O2scl is free software; you can redistribute it and/or modify
00009   it under the terms of the GNU General Public License as published by
00010   the Free Software Foundation; either version 3 of the License, or
00011   (at your option) any later version.
00012   
00013   O2scl is distributed in the hope that it will be useful,
00014   but WITHOUT ANY WARRANTY; without even the implied warranty of
00015   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016   GNU General Public License for more details.
00017   
00018   You should have received a copy of the GNU General Public License
00019   along with O2scl. If not, see <http://www.gnu.org/licenses/>.
00020 
00021   -------------------------------------------------------------------
00022 */
00023 #ifndef O2SCL_LANCZOS_H
00024 #define O2SCL_LANCZOS_H
00025 
00026 #include <cmath>
00027 #include <gsl/gsl_linalg.h>
00028 #include <o2scl/collection.h>
00029 
00030 #ifndef DOXYGENP
00031 namespace o2scl {
00032 #endif
00033 
00034   /** 
00035       \brief Lanczos diagonalization
00036       
00037       This is useful for approximating the largest eigenvalues of a
00038       symmetric matrix.
00039 
00040       The vector and matrix types can be any type which provides
00041       access via \c operator[], given suitably constructed allocation
00042       types.
00043       
00044       The tridiagonalization routine was rewritten from the EISPACK
00045       routines \c imtql1.f (but uses \c gsl_hypot() instead of \c
00046       pythag.f). 
00047 
00048       \future The function eigen_tdiag() automatically sorts the
00049       eigenvalues, which may not be necessary.
00050   */
00051   template<class vec_t, class mat_t, class alloc_vec_t, class alloc_t>
00052     class lanczos {
00053   public:
00054 
00055     lanczos() {
00056       td_iter=30;
00057       td_lasteval=0;
00058     }
00059   
00060     /** \brief Number of iterations for finding the eigenvalues of the
00061         tridiagonal matrix (default 30)
00062     */
00063     size_t td_iter;
00064 
00065     /** \brief The index for the last eigenvalue not determined if
00066         tridiagonalization fails
00067     */
00068     size_t td_lasteval;
00069 
00070     /** \brief Approximate the largest eigenvalues of a symmetric
00071         matrix \c mat using the Lanczos method
00072 
00073         Given a square matrix \c mat with size \c size, this function
00074         applies \c n_iter iterations of the Lanczos algorithm to
00075         produce \c n_iter approximate eigenvalues stored in \c
00076         eigen. As a by-product, this function also partially
00077         tridiagonalizes the matrix placing the result in \c diag and
00078         \c off_diag. Before calling this function, space must have
00079         already been allocated for \c eigen, \c diag, and \c
00080         off_diag. All three of these arrays must have at least enough
00081         space for \c n_iter elements.
00082         
00083         Choosing /c n_iter = \c size will produce all of the exact
00084         eigenvalues and the corresponding tridiagonal matrix, but this
00085         may be slower than diagonalizing the matrix directly.
00086     */
00087     int eigenvalues(size_t size, mat_t &mat, size_t n_iter, 
00088                     vec_t &eigen, vec_t &diag, vec_t &off_diag) {
00089       double t;
00090       bool cont=true;
00091       size_t i, j, k;
00092 
00093       alloc_vec_t v;
00094       alloc_vec_t w;
00095       alloc_vec_t b3;
00096       alloc_vec_t prod;
00097       alloc_t ao;
00098       ao.allocate(v,size);
00099       ao.allocate(w,size);
00100       ao.allocate(b3,size);
00101       ao.allocate(prod,size);
00102 
00103       // Pick a unit vector
00104       w[0]=1.0;
00105       for(i=1;i<size;i++) w[i]=0.0;
00106 
00107       for(i=0;i<size;i++) v[i]=0;
00108       j=0;
00109       while (cont) {
00110         if (j!=0) {
00111           for(i=0;i<size;i++) {
00112             t=w[i];
00113             w[i]=v[i]/off_diag[j-1];
00114             v[i]=-off_diag[j-1]*t;
00115           }
00116         }
00117         product(size,mat,w,prod);
00118         for(k=0;k<size;k++) v[k]+=prod[k];
00119         diag[j]=0.0;
00120         off_diag[j]=0.0;
00121 
00122         for(k=0;k<size;k++) diag[j]+=w[k]*v[k];
00123         for(k=0;k<size;k++) v[k]-=diag[j]*w[k];
00124         for(k=0;k<size;k++) off_diag[j]+=v[k]*v[k];
00125         off_diag[j]=sqrt(off_diag[j]);
00126         j++;
00127     
00128         if (j>=n_iter || off_diag[j-1]==0.0) cont=false;
00129     
00130         if (j>0) {
00131           for(k=0;k<size-1;k++) {
00132             b3[k+1]=off_diag[k];
00133             eigen[k]=diag[k];
00134           }
00135           eigen[size-1]=diag[size-1];
00136           if (eigen_tdiag(j,eigen,b3)!=0) {
00137 
00138             ao.free(v);
00139             ao.free(w);
00140             ao.free(b3);
00141             ao.free(prod);
00142 
00143             add_err_ret("Call to eigen_tdiag() in eigenvalues() failed.",
00144                         gsl_efailed);
00145           }
00146         }
00147       }
00148 
00149       ao.free(v);
00150       ao.free(w);
00151       ao.free(b3);
00152       ao.free(prod);
00153   
00154       return 0;
00155     }
00156     
00157     /** 
00158         \brief In-place diagonalization of a tri-diagonal matrix
00159 
00160         On input, the vectors \c diag and \c off_diag should both be
00161         vectors of size \c n. The diagonal entries stored in \c diag,
00162         and the \f$ n-1 \f$ off-diagonal entries should be stored in
00163         \c off_diag, starting with \c off_diag[1].  The value in \c
00164         off_diag[0] is unused. The vector \c off_diag is destroyed by
00165         the computation.
00166 
00167         This uses an implict QL method from the EISPACK routine \c
00168         imtql1. The value of \c ierr from the original Fortran routine
00169         is stored in \ref td_lasteval.
00170 
00171      */
00172     int eigen_tdiag(size_t n, vec_t &diag, vec_t &off_diag) {
00173 
00174       // 'i' is set to zero here because Cygwin complained
00175       // about uninit'ed variables. This is probably ok,
00176       // but it would be nice to double check that there
00177       // isn't a problem with setting i=0 here.
00178       int i=0,j,l,m,mml;
00179       double b,c,f,g,p,r,s,tst1,tst2;
00180 
00181       if (n==1) return 0;
00182       
00183       for(size_t ij=1;ij<n;ij++) off_diag[ij-1]=off_diag[ij];
00184       off_diag[n-1]=0.0;
00185   
00186       bool done=false;
00187   
00188       l=1;
00189       j=0;
00190   
00191       while (done==false && l<=((int)n)) {
00192     
00193         // Look for small sub-diagonal element
00194         bool idone=false;
00195         for(m=l;m<((int)n) && idone==false;m++) {
00196           tst1=fabs(diag[m-1])+fabs(diag[m]);
00197           tst2=tst1+fabs(off_diag[m-1]);
00198           if (tst2==tst1) {
00199             m--;
00200             idone=true;
00201           }
00202         }
00203     
00204         p=diag[l-1];
00205 
00206         if (m!=l && j==((int)td_iter)) {
00207           
00208           // Set error. No convergence after td_iter iterations
00209           td_lasteval=l;
00210           set_err_ret("No convergence in lanczos::eigen_tdiag()",
00211                       gsl_efailed);
00212         }
00213 
00214         if (m!=l) {
00215 
00216           j++;
00217 
00218           // Form shift
00219           g=(diag[l]-p)/(2.0*off_diag[l-1]);
00220           r=gsl_hypot(g,1.0);
00221       
00222           g=diag[m-1]-p+off_diag[l-1]/(g+(g>=0.0 ? fabs(r) : -fabs(r)));
00223           s=1.0;
00224           c=1.0;
00225           p=0.0;
00226           mml=m-l;
00227       
00228           for(int ii=1;ii<=mml;ii++) {
00229 
00230             i=m-ii;
00231             f=s*off_diag[i-1];
00232             b=c*off_diag[i-1];
00233             r=gsl_hypot(f,g);
00234             off_diag[i]=r;
00235         
00236             if (r==0.0) {
00237 
00238               // Recover from underflow
00239               diag[i]-=p;
00240               off_diag[m-1]=0.0;
00241               ii=mml+1;
00242 
00243             } else {
00244 
00245               s=f/r;
00246               c=g/r;
00247               g=diag[i]-p;
00248               r=(diag[i-1]-g)*s+2.0*c*b;
00249               p=s*r;
00250               diag[i]=g+p;
00251               g=c*r-b;
00252 
00253             }
00254 
00255           }
00256       
00257           diag[l-1]-=p;
00258           off_diag[l-1]=g;
00259           off_diag[m-1]=0.0;
00260       
00261 
00262         } else {
00263 
00264           // Order eigenvalues
00265 
00266           if (l==1) {
00267 
00268             i=1;
00269             diag[i-1]=p;
00270         
00271           } else {
00272 
00273             bool skip=false;
00274             for(int ii=2;ii<=l;ii++) {
00275               i=l+2-ii;
00276               if (p>=diag[i-2]) {
00277                 ii=l+1;
00278                 skip=true;
00279               } else {
00280                 diag[i-1]=diag[i-2];
00281               }
00282             }
00283         
00284             if (skip==false) i=1;
00285             diag[i-1]=p;
00286           }
00287 
00288           j=0;
00289           l++;
00290         }
00291     
00292       }
00293   
00294       return 0;
00295     }
00296     
00297 #ifndef DOXYGEN_INTERNAL
00298 
00299   protected:
00300       
00301     /** 
00302         \brief Naive matrix-vector product
00303         
00304         It is assumed that memory is already allocated for \c prod.
00305     */
00306     void product(size_t n, mat_t &a, vec_t &w, vec_t &prod) {
00307       size_t i, j;
00308       for(i=0;i<n;i++) {
00309         prod[i]=0.0;
00310         for(j=0;j<n;j++) {
00311           prod[i]+=a[i][j]*w[j];
00312         }
00313       }
00314       return;
00315     }
00316     
00317 #endif
00318   
00319   };
00320 
00321 #ifndef DOXYGENP
00322 }
00323 #endif
00324 
00325 #endif

Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.

Project hosting provided by SourceForge.net Logo, O2scl Sourceforge Project Page