/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: ./lib/algorithms/iterative/ImplicitlyRestartedBlockLanczos.h Copyright (C) 2015 Author: Peter Boyle Author: Chulwoo Jung Author: Yong-Chull Jang Author: Guido Cossu 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_IRBL_H #define GRID_IRBL_H #include //memset #define Glog std::cout << GridLogMessage namespace Grid { enum class LanczosType { irbl, rbl }; ///////////////////////////////////////////////////////////// // Implicitly restarted block lanczos ///////////////////////////////////////////////////////////// template class ImplicitlyRestartedBlockLanczos { private: std::string cname = std::string("ImplicitlyRestartedBlockLanczos"); int MaxIter; // Max iterations int Nstop; // Number of evecs checked for convergence int Nu; // Numbeer of vecs in the unit block int Nk; // Number of converged sought int Nm; // total number of vectors int Nblock_k; // Nk/Nu int Nblock_m; // Nm/Nu int Nconv_test_interval; // Number of skipped vectors when checking a convergence RealD eresid; IRLdiagonalisation diagonalisation; //////////////////////////////////// // Embedded objects //////////////////////////////////// SortEigen _sort; LinearOperatorBase &_Linop; OperatorFunction &_poly; ///////////////////////// // Constructor ///////////////////////// public: ImplicitlyRestartedBlockLanczos(LinearOperatorBase &Linop, // op OperatorFunction & poly, // polynomial int _Nstop, // really sought vecs int _Nconv_test_interval, // conv check interval int _Nu, // vecs in the unit block int _Nk, // sought vecs int _Nm, // total vecs RealD _eresid, // resid in lmd deficit int _MaxIter, // Max iterations IRLdiagonalisation _diagonalisation = IRLdiagonaliseWithEigen) : _Linop(Linop), _poly(poly), Nstop(_Nstop), Nconv_test_interval(_Nconv_test_interval), Nu(_Nu), Nk(_Nk), Nm(_Nm), Nblock_m(_Nm/_Nu), Nblock_k(_Nk/_Nu), //eresid(_eresid), MaxIter(10), eresid(_eresid), MaxIter(_MaxIter), diagonalisation(_diagonalisation) { assert( (Nk%Nu==0) && (Nm%Nu==0) ); }; //////////////////////////////// // Helpers //////////////////////////////// static RealD normalize(Field& v) { RealD nn = norm2(v); nn = sqrt(nn); v = v * (1.0/nn); return nn; } void orthogonalize(Field& w, std::vector& evec, int k) { typedef typename Field::scalar_type MyComplex; MyComplex ip; for(int j=0; j& evec, int k, int Nu) { typedef typename Field::scalar_type MyComplex; MyComplex ip; for(int j=0; j& eval, std::vector& evec, const std::vector& src, int& Nconv, LanczosType Impl) { switch (Impl) { case LanczosType::irbl: calc_irbl(eval,evec,src,Nconv); break; case LanczosType::rbl: calc_rbl(eval,evec,src,Nconv); break; } } void calc_irbl(std::vector& eval, std::vector& evec, const std::vector& src, int& Nconv) { std::string fname = std::string(cname+"::calc_irbl()"); GridBase *grid = evec[0]._grid; assert(grid == src[0]._grid); assert( Nu = src.size() ); Glog << std::string(74,'*') << std::endl; Glog << fname + " starting iteration 0 / "<< MaxIter<< std::endl; Glog << std::string(74,'*') << std::endl; Glog <<" -- seek Nk = "<< Nk <<" vectors"<< std::endl; Glog <<" -- accept Nstop = "<< Nstop <<" vectors"<< std::endl; Glog <<" -- total Nm = "<< Nm <<" vectors"<< std::endl; Glog <<" -- size of eval = "<< eval.size() << std::endl; Glog <<" -- size of evec = "<< evec.size() << std::endl; if ( diagonalisation == IRLdiagonaliseWithEigen ) { Glog << "Diagonalisation is Eigen "<< std::endl; } else { abort(); } Glog << std::string(74,'*') << std::endl; assert(Nm == evec.size() && Nm == eval.size()); std::vector> lmd(Nu,std::vector(Nm,0.0)); std::vector> lme(Nu,std::vector(Nm,0.0)); std::vector> lmd2(Nu,std::vector(Nm,0.0)); std::vector> lme2(Nu,std::vector(Nm,0.0)); std::vector eval2(Nm); std::vector resid(Nk); Eigen::MatrixXcd Qt = Eigen::MatrixXcd::Zero(Nm,Nm); Eigen::MatrixXcd Q = Eigen::MatrixXcd::Zero(Nm,Nm); std::vector Iconv(Nm); std::vector B(Nm,grid); // waste of space replicating std::vector f(Nu,grid); std::vector f_copy(Nu,grid); Field v(grid); Nconv = 0; RealD beta_k; // set initial vector for (int i=0; i& eval, std::vector& evec, const std::vector& src, int& Nconv) { std::string fname = std::string(cname+"::calc_rbl()"); GridBase *grid = evec[0]._grid; assert(grid == src[0]._grid); assert( Nu = src.size() ); int Np = (Nm-Nk); if (Np > 0 && MaxIter > 1) Np /= MaxIter; int Nblock_p = Np/Nu; Glog << std::string(74,'*') << std::endl; Glog << fname + " starting iteration 0 / "<< MaxIter<< std::endl; Glog << std::string(74,'*') << std::endl; Glog <<" -- seek (min) Nk = "<< Nk <<" vectors"<< std::endl; Glog <<" -- seek (inc) Np = "<< Np <<" vectors"<< std::endl; Glog <<" -- seek (max) Nm = "<< Nm <<" vectors"<< std::endl; Glog <<" -- accept Nstop = "<< Nstop <<" vectors"<< std::endl; Glog <<" -- size of eval = "<< eval.size() << std::endl; Glog <<" -- size of evec = "<< evec.size() << std::endl; if ( diagonalisation == IRLdiagonaliseWithEigen ) { Glog << "Diagonalisation is Eigen "<< std::endl; } else { abort(); } Glog << std::string(74,'*') << std::endl; assert(Nm == evec.size() && Nm == eval.size()); std::vector> lmd(Nu,std::vector(Nm,0.0)); std::vector> lme(Nu,std::vector(Nm,0.0)); std::vector> lmd2(Nu,std::vector(Nm,0.0)); std::vector> lme2(Nu,std::vector(Nm,0.0)); std::vector eval2(Nk); std::vector resid(Nm); Eigen::MatrixXcd Qt = Eigen::MatrixXcd::Zero(Nm,Nm); Eigen::MatrixXcd Q = Eigen::MatrixXcd::Zero(Nm,Nm); std::vector Iconv(Nm); std::vector B(Nm,grid); // waste of space replicating std::vector f(Nu,grid); std::vector f_copy(Nu,grid); Field v(grid); Nconv = 0; RealD beta_k; // set initial vector for (int i=0; i>& lmd, std::vector>& lme, std::vector& evec, std::vector& w, std::vector& w_copy, int b) { const RealD tiny = 1.0e-20; int Nu = w.size(); int Nm = evec.size(); assert( b < Nm/Nu ); // converts block index to full indicies for an interval [L,R) int L = Nu*b; int R = Nu*(b+1); Real beta; // 3. wk:=Avk−βkv_{k−1} for (int k=L, u=0; k0) { for (int u=0; u>& lmd, std::vector>& lme, int Nu, int Nb, int Nk, int Nm, Eigen::MatrixXcd& M) { //Glog << "packHermitBlockTriDiagMatfromEigen() begin" << '\n'; assert( Nk%Nu == 0 && Nm%Nu == 0 ); assert( Nk <= Nm ); // rearrange for ( int u=0; u QRD(Mtmp); Q = QRD.householderQ(); R = QRD.matrixQR(); // upper triangular part is the R matrix. // lower triangular part used to represent series // of Q sequence. // equivalent operation of Qprod *= Q //M = Eigen::MatrixXcd::Zero(Nm,Nm); //for (int i=0; i Nm) kmax = Nm; for (int k=i; ki) M(i,j) = conj(M(j,i)); // if (i-j > Nu || j-i > Nu) M(i,j) = 0.; // } //} //Glog << "shiftedQRDecompEigen() end" << endl; } void exampleQRDecompEigen(void) { Eigen::MatrixXd A = Eigen::MatrixXd::Zero(3,3); Eigen::MatrixXd Q = Eigen::MatrixXd::Zero(3,3); Eigen::MatrixXd R = Eigen::MatrixXd::Zero(3,3); Eigen::MatrixXd P = Eigen::MatrixXd::Zero(3,3); A(0,0) = 12.0; A(0,1) = -51.0; A(0,2) = 4.0; A(1,0) = 6.0; A(1,1) = 167.0; A(1,2) = -68.0; A(2,0) = -4.0; A(2,1) = 24.0; A(2,2) = -41.0; Glog << "matrix A before ColPivHouseholder" << std::endl; for ( int i=0; i<3; i++ ) { for ( int j=0; j<3; j++ ) { Glog << "A[" << i << "," << j << "] = " << A(i,j) << '\n'; } } Glog << std::endl; Eigen::ColPivHouseholderQR QRD(A); Glog << "matrix A after ColPivHouseholder" << std::endl; for ( int i=0; i<3; i++ ) { for ( int j=0; j<3; j++ ) { Glog << "A[" << i << "," << j << "] = " << A(i,j) << '\n'; } } Glog << std::endl; Glog << "HouseholderQ with sequence lenth = nonzeroPiviots" << std::endl; Q = QRD.householderQ().setLength(QRD.nonzeroPivots()); for ( int i=0; i<3; i++ ) { for ( int j=0; j<3; j++ ) { Glog << "Q[" << i << "," << j << "] = " << Q(i,j) << '\n'; } } Glog << std::endl; Glog << "HouseholderQ with sequence lenth = 1" << std::endl; Q = QRD.householderQ().setLength(1); for ( int i=0; i<3; i++ ) { for ( int j=0; j<3; j++ ) { Glog << "Q[" << i << "," << j << "] = " << Q(i,j) << '\n'; } } Glog << std::endl; Glog << "HouseholderQ with sequence lenth = 2" << std::endl; Q = QRD.householderQ().setLength(2); for ( int i=0; i<3; i++ ) { for ( int j=0; j<3; j++ ) { Glog << "Q[" << i << "," << j << "] = " << Q(i,j) << '\n'; } } Glog << std::endl; Glog << "matrixR" << std::endl; R = QRD.matrixR(); for ( int i=0; i<3; i++ ) { for ( int j=0; j<3; j++ ) { Glog << "R[" << i << "," << j << "] = " << R(i,j) << '\n'; } } Glog << std::endl; Glog << "rank = " << QRD.rank() << std::endl; Glog << "threshold = " << QRD.threshold() << std::endl; Glog << "matrixP" << std::endl; P = QRD.colsPermutation(); for ( int i=0; i<3; i++ ) { for ( int j=0; j<3; j++ ) { Glog << "P[" << i << "," << j << "] = " << P(i,j) << '\n'; } } Glog << std::endl; Glog << "QR decomposition without column pivoting" << std::endl; A(0,0) = 12.0; A(0,1) = -51.0; A(0,2) = 4.0; A(1,0) = 6.0; A(1,1) = 167.0; A(1,2) = -68.0; A(2,0) = -4.0; A(2,1) = 24.0; A(2,2) = -41.0; Glog << "matrix A before Householder" << std::endl; for ( int i=0; i<3; i++ ) { for ( int j=0; j<3; j++ ) { Glog << "A[" << i << "," << j << "] = " << A(i,j) << '\n'; } } Glog << std::endl; Eigen::HouseholderQR QRDplain(A); Glog << "HouseholderQ" << std::endl; Q = QRDplain.householderQ(); for ( int i=0; i<3; i++ ) { for ( int j=0; j<3; j++ ) { Glog << "Q[" << i << "," << j << "] = " << Q(i,j) << '\n'; } } Glog << std::endl; Glog << "matrix A after Householder" << std::endl; for ( int i=0; i<3; i++ ) { for ( int j=0; j<3; j++ ) { Glog << "A[" << i << "," << j << "] = " << A(i,j) << '\n'; } } Glog << std::endl; } }; } #undef Glog #endif