/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: ./lib/algorithms/iterative/ImplicitlyRestartedLanczos.h Copyright (C) 2015 Author: Peter Boyle Author: paboyle This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. See the full license in the file "LICENSE" in the top level distribution directory *************************************************************************************/ /* END LEGAL */ #ifndef GRID_IRL_H #define GRID_IRL_H #include //memset #ifdef USE_LAPACK #ifdef USE_MKL #include #else void LAPACK_dstegr(char *jobz, char *range, int *n, double *d, double *e, double *vl, double *vu, int *il, int *iu, double *abstol, int *m, double *w, double *z, int *ldz, int *isuppz, double *work, int *lwork, int *iwork, int *liwork, int *info); //#include #endif #endif #include "DenseMatrix.h" #include "EigenSort.h" namespace Grid { ///////////////////////////////////////////////////////////// // Implicitly restarted lanczos ///////////////////////////////////////////////////////////// template class ImplicitlyRestartedLanczos { const RealD small = 1.0e-16; public: int lock; int get; int Niter; int converged; int Nstop; // Number of evecs checked for convergence int Nk; // Number of converged sought int Np; // Np -- Number of spare vecs in kryloc space int Nm; // Nm -- total number of vectors RealD OrthoTime; RealD eresid; SortEigen _sort; LinearOperatorBase &_Linop; OperatorFunction &_poly; ///////////////////////// // Constructor ///////////////////////// void init(void){}; void Abort(int ff, DenseVector &evals, DenseVector > &evecs); ImplicitlyRestartedLanczos( LinearOperatorBase &Linop, // op OperatorFunction & poly, // polynmial int _Nstop, // sought vecs int _Nk, // sought vecs int _Nm, // spare vecs RealD _eresid, // resid in lmdue deficit int _Niter) : // Max iterations _Linop(Linop), _poly(poly), Nstop(_Nstop), Nk(_Nk), Nm(_Nm), eresid(_eresid), Niter(_Niter) { Np = Nm-Nk; assert(Np>0); }; ImplicitlyRestartedLanczos( LinearOperatorBase &Linop, // op OperatorFunction & poly, // polynmial int _Nk, // sought vecs int _Nm, // spare vecs RealD _eresid, // resid in lmdue deficit int _Niter) : // Max iterations _Linop(Linop), _poly(poly), Nstop(_Nk), Nk(_Nk), Nm(_Nm), eresid(_eresid), Niter(_Niter) { Np = Nm-Nk; assert(Np>0); }; ///////////////////////// // Sanity checked this routine (step) against Saad. ///////////////////////// void RitzMatrix(DenseVector& evec,int k){ if(1) return; GridBase *grid = evec[0]._grid; Field w(grid); std::cout<1 ) { if (abs(in) >1.0e-9 ) { std::cout<& lmd, DenseVector& lme, DenseVector& evec, Field& w,int Nm,int k) { assert( k< Nm ); _poly(_Linop,evec[k],w); // 3. wk:=Avk−βkv_{k−1} if(k>0){ w -= lme[k-1] * evec[k-1]; } ComplexD zalph = innerProduct(evec[k],w); // 4. αk:=(wk,vk) RealD alph = real(zalph); w = w - alph * evec[k];// 5. wk:=wk−αkvk RealD beta = normalise(w); // 6. βk+1 := ∥wk∥2. If βk+1 = 0 then Stop // 7. vk+1 := wk/βk+1 std::cout<0) { orthogonalize(w,evec,k); // orthonormalise } if(k < Nm-1) evec[k+1] = w; } void qr_decomp(DenseVector& lmd, DenseVector& lme, int Nk, int Nm, DenseVector& Qt, RealD Dsh, int kmin, int kmax) { int k = kmin-1; RealD x; RealD Fden = 1.0/hypot(lmd[k]-Dsh,lme[k]); RealD c = ( lmd[k] -Dsh) *Fden; RealD s = -lme[k] *Fden; RealD tmpa1 = lmd[k]; RealD tmpa2 = lmd[k+1]; RealD tmpb = lme[k]; lmd[k] = c*c*tmpa1 +s*s*tmpa2 -2.0*c*s*tmpb; lmd[k+1] = s*s*tmpa1 +c*c*tmpa2 +2.0*c*s*tmpb; lme[k] = c*s*(tmpa1-tmpa2) +(c*c-s*s)*tmpb; x =-s*lme[k+1]; lme[k+1] = c*lme[k+1]; for(int i=0; i& lmd, DenseVector& lme, int N1, int N2, DenseVector& Qt, GridBase *grid){ const int size = Nm; // tevals.resize(size); // tevecs.resize(size); LAPACK_INT NN = N1; double evals_tmp[NN]; double evec_tmp[NN][NN]; memset(evec_tmp[0],0,sizeof(double)*NN*NN); // double AA[NN][NN]; double DD[NN]; double EE[NN]; for (int i = 0; i< NN; i++) for (int j = i - 1; j <= i + 1; j++) if ( j < NN && j >= 0 ) { if (i==j) DD[i] = lmd[i]; if (i==j) evals_tmp[i] = lmd[i]; if (j==(i-1)) EE[j] = lme[j]; } LAPACK_INT evals_found; LAPACK_INT lwork = ( (18*NN) > (1+4*NN+NN*NN)? (18*NN):(1+4*NN+NN*NN)) ; LAPACK_INT liwork = 3+NN*10 ; LAPACK_INT iwork[liwork]; double work[lwork]; LAPACK_INT isuppz[2*NN]; char jobz = 'V'; // calculate evals & evecs char range = 'I'; // calculate all evals // char range = 'A'; // calculate all evals char uplo = 'U'; // refer to upper half of original matrix char compz = 'I'; // Compute eigenvectors of tridiagonal matrix int ifail[NN]; long long info; // int total = QMP_get_number_of_nodes(); // int node = QMP_get_node_number(); // GridBase *grid = evec[0]._grid; int total = grid->_Nprocessors; int node = grid->_processor; int interval = (NN/total)+1; double vl = 0.0, vu = 0.0; LAPACK_INT il = interval*node+1 , iu = interval*(node+1); if (iu > NN) iu=NN; double tol = 0.0; if (1) { memset(evals_tmp,0,sizeof(double)*NN); if ( il <= NN){ printf("total=%d node=%d il=%d iu=%d\n",total,node,il,iu); #ifdef USE_MKL dstegr(&jobz, &range, &NN, #else LAPACK_dstegr(&jobz, &range, &NN, #endif (double*)DD, (double*)EE, &vl, &vu, &il, &iu, // these four are ignored if second parameteris 'A' &tol, // tolerance &evals_found, evals_tmp, (double*)evec_tmp, &NN, isuppz, work, &lwork, iwork, &liwork, &info); for (int i = iu-1; i>= il-1; i--){ printf("node=%d evals_found=%d evals_tmp[%d] = %g\n",node,evals_found, i - (il-1),evals_tmp[i - (il-1)]); evals_tmp[i] = evals_tmp[i - (il-1)]; if (il>1) evals_tmp[i-(il-1)]=0.; for (int j = 0; j< NN; j++){ evec_tmp[i][j] = evec_tmp[i - (il-1)][j]; if (il>1) evec_tmp[i-(il-1)][j]=0.; } } } { // QMP_sum_double_array(evals_tmp,NN); // QMP_sum_double_array((double *)evec_tmp,NN*NN); grid->GlobalSumVector(evals_tmp,NN); grid->GlobalSumVector((double*)evec_tmp,NN*NN); } } // cheating a bit. It is better to sort instead of just reversing it, but the document of the routine says evals are sorted in increasing order. qr gives evals in decreasing order. for(int i=0;i& lmd, DenseVector& lme, int N2, int N1, DenseVector& Qt, GridBase *grid) { #ifdef USE_LAPACK const int check_lapack=0; // just use lapack if 0, check against lapack if 1 if(!check_lapack) return diagonalize_lapack(lmd,lme,N2,N1,Qt,grid); DenseVector lmd2(N1); DenseVector lme2(N1); DenseVector Qt2(N1*N1); for(int k=0; k= kmin; --j){ RealD dds = fabs(lmd[j-1])+fabs(lmd[j]); if(fabs(lme[j-1])+dds > dds){ kmax = j+1; goto continued; } } Niter = iter; #ifdef USE_LAPACK if(check_lapack){ const double SMALL=1e-8; diagonalize_lapack(lmd2,lme2,N2,N1,Qt2,grid); DenseVector lmd3(N2); for(int k=0; kSMALL) std::cout<SMALL) std::cout<SMALL) std::cout< dds){ kmin = j+1; break; } } } std::cout<& evec, int k) { double t0=-usecond()/1e6; typedef typename Field::scalar_type MyComplex; MyComplex ip; if ( 0 ) { for(int j=0; j &Qt) { for(int i=0; i K P = M − K † Compute the factorization AVM = VM HM + fM eM repeat Q=I for i = 1,...,P do QiRi =HM −θiI Q = QQi H M = Q †i H M Q i end for βK =HM(K+1,K) σK =Q(M,K) r=vK+1βK +rσK VK =VM(1:M)Q(1:M,1:K) HK =HM(1:K,1:K) →AVK =VKHK +fKe†K † Extend to an M = K + P step factorization AVM = VMHM + fMeM until convergence */ void calc(DenseVector& eval, DenseVector& evec, const Field& src, int& Nconv) { GridBase *grid = evec[0]._grid; assert(grid == src._grid); std::cout<oSites();ss++){ for(int jj=k1-1; jjoSites();ss++){ for(int jj=0; jj=Nstop ){ goto converged; } } // end of iter loop std::cout< & bq, Field &bf, DenseMatrix &H){ GridBase *grid = bq[0]._grid; RealD beta; RealD sqbt; RealD alpha; for(int i=start;i 1) std::cout< evals, DenseVector evecs){ int N= evals.size(); _sort.push(evals,evecs, evals.size(),N); } void ImplicitRestart(int TM, DenseVector &evals, DenseVector > &evecs, DenseVector &bq, Field &bf, int cont) { std::cout< H; Resize(H,Nm,Nm); #ifndef USE_LAPACK EigenSort(evals, evecs); #endif ///Assign shifts int K=Nk; int M=Nm; int P=Np; int converged=0; if(K - converged < 4) P = (M - K-1); //one // DenseVector shifts(P + shift_extra.size()); DenseVector shifts(P); for(int k = 0; k < P; ++k) shifts[k] = evals[k]; /// Shift to form a new H and q DenseMatrix Q; Resize(Q,TM,TM); Unity(Q); Shift(Q, shifts); // H is implicitly passed in in Rudy's Shift routine int ff = K; /// Shifted H defines a new K step Arnoldi factorization RealD beta = H[ff][ff-1]; RealD sig = Q[TM - 1][ff - 1]; std::cout< q Q times_real(bq, Q, TM); std::cout< &bq, Field &bf, DenseVector > & evecs,DenseVector &evals) { init(); int M=Nm; DenseMatrix H; Resize(H,Nm,Nm); Resize(evals,Nm); Resize(evecs,Nm); int ff = Lanczos_Factor(0, M, cont, bq,bf,H); // 0--M to begin with if(ff < M) { std::cout< " << it << std::endl; int lock_num = lock ? converged : 0; DenseVector tevals(M - lock_num ); DenseMatrix tevecs; Resize(tevecs,M - lock_num,M - lock_num); //check residual of polynominal TestConv(H,M, tevals, tevecs); if(converged >= Nk) break; ImplicitRestart(ff, tevals,tevecs,H); } Wilkinson(H, evals, evecs, small); // Check(); std::cout< & H,DenseMatrix &Q, DenseVector shifts) { int P; Size(shifts,P); int M; SizeSquare(Q,M); Unity(Q); int lock_num = lock ? converged : 0; RealD t_Househoulder_vector(0.0); RealD t_Househoulder_mult(0.0); for(int i=0;i ck(3), v(3); x = H[lock_num+0][lock_num+0]-shifts[i]; y = H[lock_num+1][lock_num+0]; ck[0] = x; ck[1] = y; ck[2] = 0; normalise(ck); ///Normalization cancels in PHP anyway RealD beta; Householder_vector(ck, 0, 2, v, beta); Householder_mult(H,v,beta,0,lock_num+0,lock_num+2,0); Householder_mult(H,v,beta,0,lock_num+0,lock_num+2,1); ///Accumulate eigenvector Householder_mult(Q,v,beta,0,lock_num+0,lock_num+2,1); int sw = 0; for(int k=lock_num+0;k(ck, 0, 2-sw, v, beta); Householder_mult(H,v, beta,0,k+1,k+3-sw,0); Householder_mult(H,v, beta,0,k+1,k+3-sw,1); ///Accumulate eigenvector Householder_mult(Q,v, beta,0,k+1,k+3-sw,1); } } } void TestConv(DenseMatrix & H,int SS, DenseVector &bq, Field &bf, DenseVector &tevals, DenseVector > &tevecs, int lock, int converged) { std::cout< AH; Resize(AH,SS - lock_num,SS - lock_num ); AH = GetSubMtx(H,lock_num, SS, lock_num, SS); int NN=tevals.size(); int AHsize=SS-lock_num; RealD small=1.0e-16; Wilkinson(AH, tevals, tevecs, small); #ifndef USE_LAPACK EigenSort(tevals, tevecs); #endif RealD resid_nrm= norm2(bf); if(!lock) converged = 0; #if 0 for(int i = SS - lock_num - 1; i >= SS - Nk && i >= 0; --i){ RealD diff = 0; diff = abs( tevecs[i][Nm - 1 - lock_num] ) * resid_nrm; std::cout< Q; Resize(Q,M,M); bool herm = true; Lock(H, Q, tevals[i], converged, small, SS, herm); times_real(bq, Q, bq.size()); bf = Q[M - 1][M - 1]* bf; lock_num++; } converged++; std::cout< &evals, DenseVector > &evecs) { DenseVector goodval(this->get); #ifndef USE_LAPACK EigenSort(evals,evecs); #endif int NM = Nm; DenseVector< DenseVector > V; Size(V,NM); DenseVector QZ(NM*NM); for(int i = 0; i < NM; i++){ for(int j = 0; j < NM; j++){ // evecs[i][j]; } } } /** There is some matrix Q such that for any vector y Q.e_1 = y and Q is unitary. **/ template static T orthQ(DenseMatrix &Q, DenseVector y){ int N = y.size(); //Matrix Size Fill(Q,0.0); T tau; for(int i=0;i 0.0){ T gam = conj( (y[j]/tau)/tau0 ); for(int k=0;k<=j-1;k++){ Q[k][j]=-gam*y[k]; } Q[j][j]=tau0/tau; } else { Q[j-1][j]=1.0; } tau0 = tau; } return tau; } /** There is some matrix Q such that for any vector y Q.e_k = y and Q is unitary. **/ template< class T> static T orthU(DenseMatrix &Q, DenseVector y){ T tau = orthQ(Q,y); SL(Q); return tau; } /** Wind up with a matrix with the first con rows untouched say con = 2 Q is such that Qdag H Q has {x, x, val, 0, 0, 0, 0, ...} as 1st colum and the matrix is upper hessenberg and with f and Q appropriately modidied with Q is the arnoldi factorization **/ template static void Lock(DenseMatrix &H, ///Hess mtx DenseMatrix &Q, ///Lock Transform T val, ///value to be locked int con, ///number already locked RealD small, int dfg, bool herm) { //ForceTridiagonal(H); int M = H.dim; DenseVector vec; Resize(vec,M-con); DenseMatrix AH; Resize(AH,M-con,M-con); AH = GetSubMtx(H,con, M, con, M); DenseMatrix QQ; Resize(QQ,M-con,M-con); Unity(Q); Unity(QQ); DenseVector evals; Resize(evals,M-con); DenseMatrix evecs; Resize(evecs,M-con,M-con); Wilkinson(AH, evals, evecs, small); int k=0; RealD cold = abs( val - evals[k]); for(int i=1;icon+2; j--){ DenseMatrix U; Resize(U,j-1-con,j-1-con); DenseVector z; Resize(z,j-1-con); T nm = norm(z); for(int k = con+0;k Hb; Resize(Hb,j-1-con,M); for(int a = 0;a Qb; Resize(Qb,M,M); for(int a = 0;a Hc; Resize(Hc,M,M); for(int a = 0;a