36#define Glog std::cout << GridLogMessage
46 std::string
cname = std::string(
"ImplicitlyRestartedBlockLanczosCoarse");
71 void VectorPoly(std::vector<Field> &in,std::vector<Field> &out)
75 for(
int r=0;r<in.size();r+=
mrhs){
76 for(
int rr=0;rr<
mrhs;rr++){
78 if(rrr >= in.size()) rrr = 0;
82 for(
int rr=0;rr<
mrhs;rr++){
94 for(
int rr=0;rr<
mrhs;rr++){
97 _Linop.HermOp(mrhs_in,mrhs_out);
110 int _Nconv_test_interval,
116 IRBLdiagonalisation _diagonalisation = IRBLdiagonaliseWithEigen)
119 Nu(_Nu),
Nk(_Nk),
Nm(_Nm),
124 { assert( (
Nk%
Nu==0) && (
Nm%
Nu==0) ); };
137 void orthogonalize(Field& w, std::vector<Field>& evec,
int k,
int if_print=0)
139 typedef typename Field::scalar_type MyComplex;
142 for(
int j=0; j<k; ++j){
145 if( norm(ip)/
norm2(w) > 1e-14)
146 Glog<<
"orthogonalize before: "<<j<<
" of "<<k<<
" "<< ip <<std::endl;
147 w = w - ip * evec[j];
150 if( norm(ip)/
norm2(w) > 1e-14)
151 Glog<<
"orthogonalize after: "<<j<<
" of "<<k<<
" "<< ip <<std::endl;
161 void orthogonalize(std::vector<Field>& w,
int _Nu, std::vector<Field>& evec,
int k,
int if_print=0)
163 typedef typename Field::scalar_type MyComplex;
166 for(
int j=0; j<k; ++j){
167 for(
int i=0; i<_Nu; ++i){
169 w[i] = w[i] - ip * evec[j];
171 for(
int i=0; i<_Nu; ++i)
177 typedef typename Field::scalar_type MyComplex;
180 for(
int j=0; j<k; ++j){
182 w = w - ip * evec[j*
Nu];
188 std::vector<Field>& evec,
189 const std::vector<Field>& src,
int& Nconv, LanczosType Impl)
192 case LanczosType::irbl:
196 case LanczosType::rbl:
203 std::vector<Field>& evec,
204 const std::vector<Field>& src,
int& Nconv)
206 std::string fname = std::string(
cname+
"::calc_irbl()");
208 assert(grid == src[0].
Grid());
209 assert(
Nu = src.size() );
211 Glog << std::string(74,
'*') << std::endl;
212 Glog << fname +
" starting iteration 0 / "<<
MaxIter<< std::endl;
213 Glog << std::string(74,
'*') << std::endl;
214 Glog <<
" -- seek Nk = "<<
Nk <<
" vectors"<< std::endl;
215 Glog <<
" -- accept Nstop = "<<
Nstop <<
" vectors"<< std::endl;
216 Glog <<
" -- total Nm = "<<
Nm <<
" vectors"<< std::endl;
217 Glog <<
" -- size of eval = "<<
eval.size() << std::endl;
218 Glog <<
" -- size of evec = "<< evec.size() << std::endl;
220 Glog <<
"Diagonalisation is Eigen "<< std::endl;
223 Glog <<
"Diagonalisation is LAPACK "<< std::endl;
228 Glog << std::string(74,
'*') << std::endl;
230 assert(
Nm == evec.size() &&
Nm ==
eval.size());
232 std::vector<std::vector<ComplexD>> lmd(
Nu,std::vector<ComplexD>(
Nm,0.0));
233 std::vector<std::vector<ComplexD>> lme(
Nu,std::vector<ComplexD>(
Nm,0.0));
234 std::vector<std::vector<ComplexD>> lmd2(
Nu,std::vector<ComplexD>(
Nm,0.0));
235 std::vector<std::vector<ComplexD>> lme2(
Nu,std::vector<ComplexD>(
Nm,0.0));
236 std::vector<RealD> eval2(
Nm);
237 std::vector<RealD> resid(
Nk);
239 Eigen::MatrixXcd Qt = Eigen::MatrixXcd::Zero(
Nm,
Nm);
240 Eigen::MatrixXcd Q = Eigen::MatrixXcd::Zero(
Nm,
Nm);
242 std::vector<int> Iconv(
Nm);
243 std::vector<Field>
B(
Nm,grid);
245 std::vector<Field> f(
Nu,grid);
246 std::vector<Field> f_copy(
Nu,grid);
254 for (
int i=0; i<
Nu; ++i) {
255 Glog <<
"norm2(src[" << i <<
"])= "<<
norm2(src[i]) << std::endl;
266 for(iter = 0; iter<
MaxIter; ++iter){
268 Glog <<
"#Restart iteration = "<< iter << std::endl;
273 for(
int u=0; u<
Nu; ++u){
274 for(
int k=0; k<
Nm; ++k){
275 lmd2[u][k] = lmd[u][k];
276 lme2[u][k] = lme[u][k];
279 Qt = Eigen::MatrixXcd::Identity(
Nm,
Nm);
283 for(
int i=0; i<
Nm; ++i){
295 Eigen::MatrixXcd BTDM = Eigen::MatrixXcd::Identity(
Nm,
Nm);
296 Q = Eigen::MatrixXcd::Identity(
Nm,
Nm);
300 for(
int ip=
Nk; ip<
Nm; ++ip){
301 Glog <<
" ip "<<ip<<
" / "<<
Nm<<std::endl;
307 for(
int i=0; i<k2; ++i)
B[i] = 0.0;
308 for(
int j=0; j<k2; ++j){
309 for(
int k=0; k<
Nm; ++k){
310 B[j].Checkerboard() = evec[k].Checkerboard();
311 B[j] += evec[k]*Q(k,j);
314 for(
int i=0; i<k2; ++i) evec[i] =
B[i];
320 for(
int u=0; u<
Nu; ++u){
321 for(
int k=0; k<
Nm; ++k){
322 lmd2[u][k] = lmd[u][k];
323 lme2[u][k] = lme[u][k];
326 Qt = Eigen::MatrixXcd::Identity(
Nm,
Nm);
330 for(
int i=0; i<
Nk; ++i){
339 Glog <<
" #Convergence test: "<<std::endl;
340 for(
int k = 0; k<
Nk; ++k)
B[k]=0.0;
341 for(
int j = 0; j<
Nk; ++j){
342 for(
int k = 0; k<
Nk; ++k){
343 B[j].Checkerboard() = evec[k].Checkerboard();
344 B[j] += evec[k]*Qt(k,j);
349 for(
int i=0; i<
Nk; ++i){
356 eval2[i] = vnum/vden;
361 std::cout.precision(13);
362 std::cout <<
"[" << std::setw(4)<< std::setiosflags(std::ios_base::right) <<i<<
"] ";
363 std::cout <<
"eval = "<<std::setw(20)<< std::setiosflags(std::ios_base::left)<< eval2[i];
364 std::cout <<
" resid^2 = "<< std::setw(20)<< std::setiosflags(std::ios_base::right)<< vv<< std::endl;
375 Glog <<
" #modes converged: "<<Nconv<<std::endl;
376 for(
int i=0; i<Nconv; ++i){
377 std::cout.precision(13);
378 std::cout <<
"[" << std::setw(4)<< std::setiosflags(std::ios_base::right) <<Iconv[i]<<
"] ";
379 std::cout <<
"eval_conv = "<<std::setw(20)<< std::setiosflags(std::ios_base::left)<< eval2[Iconv[i]];
380 std::cout <<
" resid^2 = "<< std::setw(20)<< std::setiosflags(std::ios_base::right)<< resid[Iconv[i]]<< std::endl;
383 if ( Nconv>=
Nstop )
break;
387 Glog << std::string(74,
'*') << std::endl;
389 Glog << fname +
" NOT converged ; Summary :\n";
391 Glog << fname +
" CONVERGED ; Summary :\n";
394 evec.resize(Nconv,grid);
395 for(
int i=0; i<Nconv; ++i){
396 eval[i] = eval2[Iconv[i]];
397 evec[i] =
B[Iconv[i]];
401 Glog << std::string(74,
'*') << std::endl;
402 Glog <<
" -- Iterations = "<< iter <<
"\n";
403 Glog <<
" -- beta(k) = "<< beta_k <<
"\n";
404 Glog <<
" -- Nconv = "<< Nconv <<
"\n";
405 Glog << std::string(74,
'*') << std::endl;
411 std::vector<Field>& evec,
412 const std::vector<Field>& src,
int& Nconv)
414 std::string fname = std::string(
cname+
"::calc_rbl()");
416 assert(grid == src[0].
Grid());
417 assert(
Nu = src.size() );
421 int Nblock_p = Np/
Nu;
424 Glog << std::string(74,
'*') << std::endl;
425 Glog << fname +
" starting iteration 0 / "<<
MaxIter<< std::endl;
426 Glog << std::string(74,
'*') << std::endl;
427 Glog <<
" -- seek (min) Nk = "<<
Nk <<
" vectors"<< std::endl;
428 Glog <<
" -- seek (inc) Np = "<< Np <<
" vectors"<< std::endl;
429 Glog <<
" -- seek (max) Nm = "<<
Nm <<
" vectors"<< std::endl;
430 Glog <<
" -- accept Nstop = "<<
Nstop <<
" vectors"<< std::endl;
431 Glog <<
" -- size of eval = "<<
eval.size() << std::endl;
432 Glog <<
" -- size of evec = "<< evec.size() << std::endl;
434 Glog <<
"Diagonalisation is Eigen "<< std::endl;
437 Glog <<
"Diagonalisation is LAPACK "<< std::endl;
442 Glog << std::string(74,
'*') << std::endl;
444 assert(
Nm == evec.size() &&
Nm ==
eval.size());
446 std::vector<std::vector<ComplexD>> lmd(
Nu,std::vector<ComplexD>(
Nm,0.0));
447 std::vector<std::vector<ComplexD>> lme(
Nu,std::vector<ComplexD>(
Nm,0.0));
448 std::vector<std::vector<ComplexD>> lmd2(
Nu,std::vector<ComplexD>(
Nm,0.0));
449 std::vector<std::vector<ComplexD>> lme2(
Nu,std::vector<ComplexD>(
Nm,0.0));
450 std::vector<RealD> eval2(
Nk);
451 std::vector<RealD> resid(
Nm);
453 Eigen::MatrixXcd Qt = Eigen::MatrixXcd::Zero(
Nm,
Nm);
454 Eigen::MatrixXcd Q = Eigen::MatrixXcd::Zero(
Nm,
Nm);
456 std::vector<int> Iconv(
Nm);
459 std::vector<Field>
B(1,grid);
461 std::vector<Field> f(
Nu,grid);
462 std::vector<Field> f_copy(
Nu,grid);
470 for (
int i=0; i<
Nu; ++i) {
471 Glog <<
"norm2(src[" << i <<
"])= "<<
norm2(src[i]) << std::endl;
474 Glog <<
"norm2(evec[" << i <<
"])= "<<
norm2(evec[i]) << std::endl;
483 int Nblock_l, Nblock_r;
487 for(iter = 0; iter<
MaxIter; ++iter){
489 Glog <<
"#Restart iteration = "<< iter << std::endl;
491 Nblock_l =
Nblock_k + iter*Nblock_p;
492 Nblock_r = Nblock_l + Nblock_p;
498 for(
int b=Nblock_l; b<Nblock_r; ++b)
blockwiseStep(lmd,lme,evec,f,f_copy,b);
501 for(
int u=0; u<
Nu; ++u){
502 for(
int k=0; k<Nr; ++k){
503 lmd2[u][k] = lmd[u][k];
504 lme2[u][k] = lme[u][k];
507 Qt = Eigen::MatrixXcd::Identity(Nr,Nr);
509 _sort.push(eval2,Nr);
510 Glog <<
"#Ritz value: "<< std::endl;
511 for(
int i=0; i<Nr; ++i){
512 std::cout.precision(13);
513 std::cout <<
"[" << std::setw(4)<< std::setiosflags(std::ios_base::right) <<i<<
"] ";
514 std::cout <<
"Rval = "<<std::setw(20)<< std::setiosflags(std::ios_base::left)<< eval2[i] << std::endl;
518 Glog <<
" #Convergence test: "<<std::endl;
523 Glog <<
" #rotation for next check point evec"
524 << std::setw(4)<< std::setiosflags(std::ios_base::right)
525 <<
"["<< j <<
"]" <<std::endl;
526 for(
int k = 0; k<Nr; ++k){
527 B[0].Checkerboard() = evec[k].Checkerboard();
528 B[0] += evec[k]*Qt(k,j);
535 eval2[j] = vnum/vden;
540 std::cout.precision(13);
541 std::cout <<
"[" << std::setw(4)<< std::setiosflags(std::ios_base::right) <<j<<
"] ";
542 std::cout <<
"eval = "<<std::setw(20)<< std::setiosflags(std::ios_base::left)<< eval2[j];
543 std::cout <<
" resid^2 = "<< std::setw(20)<< std::setiosflags(std::ios_base::right)<< vv<< std::endl;
556 Glog <<
" #modes converged: "<<Nconv<<std::endl;
557 for(
int i=0; i<Nconv; ++i){
558 std::cout.precision(13);
559 std::cout <<
"[" << std::setw(4)<< std::setiosflags(std::ios_base::right) <<Iconv[i]<<
"] ";
560 std::cout <<
"eval_conv = "<<std::setw(20)<< std::setiosflags(std::ios_base::left)<< eval2[Iconv[i]];
561 std::cout <<
" resid^2 = "<< std::setw(20)<< std::setiosflags(std::ios_base::right)<< resid[Iconv[i]]<< std::endl;
565 if ( Nconv_guess >=
Nstop )
break;
569 Glog << std::string(74,
'*') << std::endl;
570 if ( Nconv_guess <
Nstop ) {
571 Glog << fname +
" NOT converged ; Summary :\n";
573 Glog << fname +
" CONVERGED ; Summary :\n";
576 std::vector<Field> Btmp(
Nstop,grid);
578 for(
int i=0; i<
Nstop; ++i){
580 for(
int k = 0; k<Nr; ++k){
581 Btmp[i].Checkerboard() = evec[k].Checkerboard();
582 Btmp[i] += evec[k]*Qt(k,i);
589 v -= vnum/vden*Btmp[i];
593 std::cout.precision(13);
594 std::cout <<
"[" << std::setw(4)<< std::setiosflags(std::ios_base::right) <<i<<
"] ";
595 std::cout <<
"eval = "<<std::setw(20)<< std::setiosflags(std::ios_base::left)<< vnum/vden;
596 std::cout <<
" resid^2 = "<< std::setw(20)<< std::setiosflags(std::ios_base::right)<< vv<< std::endl;
599 for(
int i=0; i<
Nstop; ++i) evec[i] = Btmp[i];
601 evec.resize(
Nstop,grid);
604 Glog << std::string(74,
'*') << std::endl;
605 Glog <<
" -- Iterations = "<< iter <<
"\n";
607 Glog <<
" -- Nconv = "<< Nconv <<
"\n";
608 Glog <<
" -- Nconv (guess) = "<< Nconv_guess <<
"\n";
609 Glog << std::string(74,
'*') << std::endl;
615 std::vector<std::vector<ComplexD>>& lme,
616 std::vector<Field>& evec,
617 std::vector<Field>& w,
618 std::vector<Field>& w_copy,
621 const RealD tiny = 1.0e-20;
624 int Nm = evec.size();
639 while ( k_start <
Nu) {
640 Glog <<
"k_start= "<<k_start<<
" Poly["<<L+k_start<<
","<<L+k_start+
mrhs-1<<
"]"<<std::endl;
641 for (
int u=0; u<
mrhs; ++u) in[u] = evec[L+k_start+u];
643 for (
int u=0; u<
mrhs; ++u) w[k_start+u] = out[u];
650 for (
int u=0; u<
Nu; ++u) {
652 for (
int k=L-
Nu+u; k<L; ++k) {
653 w[u] = w[u] - evec[k] *
conjugate(lme[u][k]);
665 for (
int u=0; u<
Nu; ++u) {
666 for (
int k=L+u; k<R; ++k) {
671 lmd[u][L+u] =
real(lmd[u][L+u]);
675 for (
int u=0; u<
Nu; ++u) {
676 for (
int k=L; k<R; ++k) {
677 w[u] = w[u] - evec[k]*lmd[u][k];
687 for (
int u=0; u<
Nu; ++u) {
688 for (
int k=L; k<R; ++k) {
697 for (
int u=1; u<
Nu; ++u) {
703 for (
int u=0; u<
Nu; ++u) {
705 for (
int v=u; v<
Nu; ++v) {
708 lme[u][L+u] =
real(lme[u][L+u]);
712 for (
int u=0; u<
Nu; ++u) {
714 assert (!isnan(
norm2(w[u])));
715 for (
int k=L+u; k<R; ++k) {
722 for (
int u=0; u<
Nu; ++u) {
731 std::vector<std::vector<ComplexD>>& lmd,
732 std::vector<std::vector<ComplexD>>& lme,
734 Eigen::MatrixXcd & Qt,
739 Eigen::MatrixXcd BlockTriDiag = Eigen::MatrixXcd::Zero(
Nk,
Nk);
741 for (
int u=0; u<
Nu; ++u ) {
742 for (
int k=0; k<
Nk; ++k ) {
743 BlockTriDiag(k,u+(k/
Nu)*
Nu) = lmd[u][k];
747 for (
int u=0; u<
Nu; ++u ) {
748 for (
int k=
Nu; k<
Nk; ++k ) {
750 BlockTriDiag(u+(k/
Nu)*
Nu,k-
Nu) = lme[u][k-
Nu];
755 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXcd> eigensolver(BlockTriDiag);
757 for (
int i = 0; i <
Nk; i++) {
758 eval[
Nk-1-i] = eigensolver.eigenvalues()(i);
760 for (
int i = 0; i <
Nk; i++) {
761 for (
int j = 0; j <
Nk; j++) {
762 Qt(j,
Nk-1-i) = eigensolver.eigenvectors()(j,i);
770 void diagonalize_lapack(std::vector<RealD>&
eval,
771 std::vector<std::vector<ComplexD>>& lmd,
772 std::vector<std::vector<ComplexD>>& lme,
774 Eigen::MatrixXcd & Qt,
777 Glog <<
"diagonalize_lapack: Nu= "<<
Nu<<
" Nk= "<<
Nk<<
" Nm= "<<std::endl;
780 Eigen::MatrixXcd BlockTriDiag = Eigen::MatrixXcd::Zero(
Nk,
Nk);
782 for (
int u=0; u<
Nu; ++u ) {
783 for (
int k=0; k<
Nk; ++k ) {
785 BlockTriDiag(k,u+(k/
Nu)*
Nu) = lmd[u][k];
789 for (
int u=0; u<
Nu; ++u ) {
790 for (
int k=
Nu; k<
Nk; ++k ) {
793 BlockTriDiag(u+(k/
Nu)*
Nu,k-
Nu) = lme[u][k-
Nu];
802 double *evals_tmp = (
double *) malloc(NN*
sizeof(
double));
803 MKL_Complex16 *evec_tmp = (MKL_Complex16 *) malloc(NN*NN*
sizeof(MKL_Complex16));
804 MKL_Complex16 *DD = (MKL_Complex16 *) malloc(NN*NN*
sizeof(MKL_Complex16));
805 for (
int i = 0; i< NN; i++) {
806 for (
int j = 0; j <NN ; j++) {
807 evec_tmp[i*NN+j].real=0.;
808 evec_tmp[i*NN+j].imag=0.;
809 DD[i*NN+j].real=BlockTriDiag(i,j).real();
810 DD[i*NN+j].imag=BlockTriDiag(i,j).imag();
814 MKL_INT lwork = (3*NN);
815 MKL_INT lrwork = (24*NN);
816 MKL_INT liwork = NN*10 ;
817 MKL_INT iwork[liwork];
818 double rwork[lrwork];
820 MKL_Complex16 *work = (MKL_Complex16 *) malloc(lwork*
sizeof(MKL_Complex16));
821 MKL_INT isuppz[2*NN];
831 int interval = (NN/total)+1;
832 double vl = 0.0, vu = 0.0;
833 MKL_INT il = interval*node+1 , iu = interval*(node+1);
835 Glog <<
"node "<<node<<
"il "<<il<<
"iu "<<iu<<std::endl;
838 memset(evals_tmp,0,
sizeof(
double)*NN);
840 zheevr(&jobz, &range, &uplo, &NN,
844 &evals_found, evals_tmp, (MKL_Complex16*)evec_tmp, &NN,
851 for (
int i = iu-1; i>= il-1; i--){
852 evals_tmp[i] = evals_tmp[i - (il-1)];
853 if (il>1) evals_tmp[i-(il-1)]=0.;
854 for (
int j = 0; j< NN; j++){
855 evec_tmp[i*NN+j] = evec_tmp[(i - (il-1))*NN+j];
857 evec_tmp[(i-(il-1))*NN+j].
imag=0.;
858 evec_tmp[(i-(il-1))*NN+j].
real=0.;
868 for (
int i = 0; i <
Nk; i++)
869 eval[
Nk-1-i] = evals_tmp[i];
870 for (
int i = 0; i <
Nk; i++) {
871 for (
int j = 0; j <
Nk; j++) {
873 Qt(j,
Nk-1-i)=std::complex<double>
883 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXcd> eigensolver(BlockTriDiag);
885 for (
int i = 0; i <
Nk; i++) {
886 Glog <<
"eval = "<<i<<
" " <<
eval[
Nk-1-i] <<
" "<< eigensolver.eigenvalues()(i) <<std::endl;
889 for (
int i = 0; i <
Nk; i++) {
890 for (
int j = 0; j <
Nk; j++) {
893 MKL_Complex16 tmp = evec_tmp[i*
Nk+j];
897 Glog<<
"Qt "<<j<<
" "<<
Nk-1-i<<
" = " << norm(Qt(j,
Nk-1-i))<<
" "<<
898 norm(eigensolver.eigenvectors()(j,i)) <<std::endl;
913 std::vector<std::vector<ComplexD>>& lmd,
914 std::vector<std::vector<ComplexD>>& lme,
916 Eigen::MatrixXcd & Qt,
919 Qt = Eigen::MatrixXcd::Identity(
Nm,
Nm);
924 diagonalize_lapack(
eval,lmd,lme,
Nu,
Nk,
Nm,Qt,grid);
933 std::vector<std::vector<ComplexD>>& lmd,
934 std::vector<std::vector<ComplexD>>& lme,
935 int Nu,
int Nb,
int Nk,
int Nm,
941 M = Eigen::MatrixXcd::Zero(
Nk,
Nk);
944 for (
int u=0; u<
Nu; ++u ) {
945 for (
int k=0; k<
Nk; ++k ) {
946 M(k,u+(k/
Nu)*
Nu) = lmd[u][k];
950 for (
int u=0; u<
Nu; ++u ) {
951 for (
int k=
Nu; k<
Nk; ++k ) {
961 std::vector<std::vector<ComplexD>>& lmd,
962 std::vector<std::vector<ComplexD>>& lme,
963 int Nu,
int Nb,
int Nk,
int Nm,
971 for (
int u=0; u<
Nu; ++u ) {
972 for (
int k=0; k<
Nk; ++k ) {
973 lmd[u][k] = M(k,u+(k/
Nu)*
Nu);
977 for (
int u=0; u<
Nu; ++u ) {
978 for (
int k=
Nu; k<
Nk; ++k ) {
989 Eigen::MatrixXcd& Qprod)
992 Eigen::MatrixXcd Q = Eigen::MatrixXcd::Zero(
Nm,
Nm);
993 Eigen::MatrixXcd R = Eigen::MatrixXcd::Zero(
Nm,
Nm);
994 Eigen::MatrixXcd Mtmp = Eigen::MatrixXcd::Zero(
Nm,
Nm);
997 for (
int i=0; i<
Nm; ++i ) {
998 Mtmp(i,i) = M(i,i) - Dsh;
1001 Eigen::HouseholderQR<Eigen::MatrixXcd> QRD(Mtmp);
1002 Q = QRD.householderQ();
1026 Mtmp = Eigen::MatrixXcd::Zero(
Nm,
Nm);
1029 for (
int i=0; i<
Nm; ++i) {
1030 for (
int j=0; j<
Nm-(
Nu+1); ++j) {
1031 for (
int k=0; k<
Nu+1+j; ++k) {
1032 Mtmp(i,j) += Qprod(i,k)*Q(k,j);
1037 for (
int i=0; i<
Nm; ++i) {
1038 for (
int j=
Nm-(
Nu+1); j<
Nm; ++j) {
1039 for (
int k=0; k<
Nm; ++k) {
1040 Mtmp(i,j) += Qprod(i,k)*Q(k,j);
1057 Mtmp = Eigen::MatrixXcd::Zero(
Nm,
Nm);
1059 for (
int a=0, i=0, kmax=0; a<
Nu+1; ++a) {
1060 for (
int j=0; j<
Nm-a; ++j) {
1063 if (kmax >
Nm) kmax =
Nm;
1064 for (
int k=i; k<kmax; ++k) {
1065 Mtmp(i,j) += R(i,k)*Q(k,j);
1067 Mtmp(j,i) =
conj(Mtmp(i,j));
1072 for (
int i=0; i<
Nm; ++i) {
1073 Mtmp(i,i) =
real(Mtmp(i,i)) + Dsh;
1093 Eigen::MatrixXd A = Eigen::MatrixXd::Zero(3,3);
1094 Eigen::MatrixXd Q = Eigen::MatrixXd::Zero(3,3);
1095 Eigen::MatrixXd R = Eigen::MatrixXd::Zero(3,3);
1096 Eigen::MatrixXd P = Eigen::MatrixXd::Zero(3,3);
1108 Glog <<
"matrix A before ColPivHouseholder" << std::endl;
1109 for (
int i=0; i<3; i++ ) {
1110 for (
int j=0; j<3; j++ ) {
1111 Glog <<
"A[" << i <<
"," << j <<
"] = " << A(i,j) <<
'\n';
1116 Eigen::ColPivHouseholderQR<Eigen::MatrixXd> QRD(A);
1118 Glog <<
"matrix A after ColPivHouseholder" << std::endl;
1119 for (
int i=0; i<3; i++ ) {
1120 for (
int j=0; j<3; j++ ) {
1121 Glog <<
"A[" << i <<
"," << j <<
"] = " << A(i,j) <<
'\n';
1126 Glog <<
"HouseholderQ with sequence lenth = nonzeroPiviots" << std::endl;
1127 Q = QRD.householderQ().setLength(QRD.nonzeroPivots());
1128 for (
int i=0; i<3; i++ ) {
1129 for (
int j=0; j<3; j++ ) {
1130 Glog <<
"Q[" << i <<
"," << j <<
"] = " << Q(i,j) <<
'\n';
1135 Glog <<
"HouseholderQ with sequence lenth = 1" << std::endl;
1136 Q = QRD.householderQ().setLength(1);
1137 for (
int i=0; i<3; i++ ) {
1138 for (
int j=0; j<3; j++ ) {
1139 Glog <<
"Q[" << i <<
"," << j <<
"] = " << Q(i,j) <<
'\n';
1144 Glog <<
"HouseholderQ with sequence lenth = 2" << std::endl;
1145 Q = QRD.householderQ().setLength(2);
1146 for (
int i=0; i<3; i++ ) {
1147 for (
int j=0; j<3; j++ ) {
1148 Glog <<
"Q[" << i <<
"," << j <<
"] = " << Q(i,j) <<
'\n';
1153 Glog <<
"matrixR" << std::endl;
1155 for (
int i=0; i<3; i++ ) {
1156 for (
int j=0; j<3; j++ ) {
1157 Glog <<
"R[" << i <<
"," << j <<
"] = " << R(i,j) <<
'\n';
1162 Glog <<
"rank = " << QRD.rank() << std::endl;
1163 Glog <<
"threshold = " << QRD.threshold() << std::endl;
1165 Glog <<
"matrixP" << std::endl;
1166 P = QRD.colsPermutation();
1167 for (
int i=0; i<3; i++ ) {
1168 for (
int j=0; j<3; j++ ) {
1169 Glog <<
"P[" << i <<
"," << j <<
"] = " << P(i,j) <<
'\n';
1175 Glog <<
"QR decomposition without column pivoting" << std::endl;
1187 Glog <<
"matrix A before Householder" << std::endl;
1188 for (
int i=0; i<3; i++ ) {
1189 for (
int j=0; j<3; j++ ) {
1190 Glog <<
"A[" << i <<
"," << j <<
"] = " << A(i,j) <<
'\n';
1195 Eigen::HouseholderQR<Eigen::MatrixXd> QRDplain(A);
1197 Glog <<
"HouseholderQ" << std::endl;
1198 Q = QRDplain.householderQ();
1199 for (
int i=0; i<3; i++ ) {
1200 for (
int j=0; j<3; j++ ) {
1201 Glog <<
"Q[" << i <<
"," << j <<
"] = " << Q(i,j) <<
'\n';
1206 Glog <<
"matrix A after Householder" << std::endl;
1207 for (
int i=0; i<3; i++ ) {
1208 for (
int j=0; j<3; j++ ) {
1209 Glog <<
"A[" << i <<
"," << j <<
"] = " << A(i,j) <<
'\n';
accelerator_inline Grid_simd< S, V > sqrt(const Grid_simd< S, V > &r)
accelerator_inline sobj eval(const uint64_t ss, const sobj &arg)
Lattice< vobj > real(const Lattice< vobj > &lhs)
Lattice< vobj > imag(const Lattice< vobj > &lhs)
Lattice< vobj > conjugate(const Lattice< vobj > &lhs)
ComplexD innerProduct(const Lattice< vobj > &left, const Lattice< vobj > &right)
RealD norm2(const Lattice< vobj > &arg)
void InsertSlice(const Lattice< vobj > &lowDim, Lattice< vobj > &higherDim, int slice, int orthog)
void ExtractSlice(Lattice< vobj > &lowDim, const Lattice< vobj > &higherDim, int slice, int orthog)
#define NAMESPACE_BEGIN(A)
std::complex< RealD > ComplexD
void GlobalSumVector(RealF *, int N)
void orthogonalize(Field &w, std::vector< Field > &evec, int k, int if_print=0)
static RealD normalize(Field &v, int if_print=0)
OperatorFunction< Field > & _poly
ImplicitlyRestartedBlockLanczosCoarse(LinearOperatorBase< Field > &Linop, GridBase *f_Grid, GridBase *mrhs_Grid, int _mrhs, OperatorFunction< Field > &poly, int _Nstop, int _Nconv_test_interval, int _Nu, int _Nk, int _Nm, RealD _eresid, int _MaxIter, IRBLdiagonalisation _diagonalisation=IRBLdiagonaliseWithEigen)
void SingleOperator(Field &in, Field &out)
void diagonalize(std::vector< RealD > &eval, std::vector< std::vector< ComplexD > > &lmd, std::vector< std::vector< ComplexD > > &lme, int Nu, int Nk, int Nm, Eigen::MatrixXcd &Qt, GridBase *grid)
void blockwiseStep(std::vector< std::vector< ComplexD > > &lmd, std::vector< std::vector< ComplexD > > &lme, std::vector< Field > &evec, std::vector< Field > &w, std::vector< Field > &w_copy, int b)
void reorthogonalize(Field &w, std::vector< Field > &evec, int k)
void calc_rbl(std::vector< RealD > &eval, std::vector< Field > &evec, const std::vector< Field > &src, int &Nconv)
void orthogonalize(std::vector< Field > &w, int _Nu, std::vector< Field > &evec, int k, int if_print=0)
void diagonalize_Eigen(std::vector< RealD > &eval, std::vector< std::vector< ComplexD > > &lmd, std::vector< std::vector< ComplexD > > &lme, int Nu, int Nk, int Nm, Eigen::MatrixXcd &Qt, GridBase *grid)
void VectorPoly(std::vector< Field > &in, std::vector< Field > &out)
void shiftedQRDecompEigen(Eigen::MatrixXcd &M, int Nu, int Nm, RealD Dsh, Eigen::MatrixXcd &Qprod)
void unpackHermitBlockTriDiagMatToEigen(std::vector< std::vector< ComplexD > > &lmd, std::vector< std::vector< ComplexD > > &lme, int Nu, int Nb, int Nk, int Nm, Eigen::MatrixXcd &M)
void exampleQRDecompEigen(void)
void orthogonalize_blockhead(Field &w, std::vector< Field > &evec, int k, int Nu)
LinearOperatorBase< Field > & _Linop
void calc(std::vector< RealD > &eval, std::vector< Field > &evec, const std::vector< Field > &src, int &Nconv, LanczosType Impl)
IRBLdiagonalisation diagonalisation
void calc_irbl(std::vector< RealD > &eval, std::vector< Field > &evec, const std::vector< Field > &src, int &Nconv)
void packHermitBlockTriDiagMatfromEigen(std::vector< std::vector< ComplexD > > &lmd, std::vector< std::vector< ComplexD > > &lme, int Nu, int Nb, int Nk, int Nm, Eigen::MatrixXcd &M)