Grid 0.7.0
ImplicitlyRestartedBlockLanczosCoarse.h
Go to the documentation of this file.
1 /*************************************************************************************
2
3 Grid physics library, www.github.com/paboyle/Grid
4
5 Source file: ./lib/algorithms/iterative/ImplicitlyRestartedBlockLanczosCoarse.h
6
7 Copyright (C) 2023
8
9Author: Peter Boyle <paboyle@ph.ed.ac.uk>
10Author: Yong-Chull Jang <ypj@quark.phy.bnl.gov>
11Author: Chulwoo Jung <chulwoo@bnl.gov>
12
13 This program is free software; you can redistribute it and/or modify
14 it under the terms of the GNU General Public License as published by
15 the Free Software Foundation; either version 2 of the License, or
16 (at your option) any later version.
17
18 This program is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU General Public License for more details.
22
23 You should have received a copy of the GNU General Public License along
24 with this program; if not, write to the Free Software Foundation, Inc.,
25 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26
27 See the full license in the file "LICENSE" in the top level distribution directory
28 *************************************************************************************/
29 /* END LEGAL */
30#pragma once
31
32#include <string.h> //memset
33
35
36#define Glog std::cout << GridLogMessage
37
39// Implicitly restarted block lanczos
41template<class Field>
43
44private:
45
46 std::string cname = std::string("ImplicitlyRestartedBlockLanczosCoarse");
47 int MaxIter; // Max iterations
48 int Nstop; // Number of evecs checked for convergence
49 int Nu; // Number of vecs in the unit block
50 int Nk; // Number of converged sought
51 int Nm; // total number of vectors
52 int Nblock_k; // Nk/Nu
53 int Nblock_m; // Nm/Nu
54 int Nconv_test_interval; // Number of skipped vectors when checking a convergence
56 IRBLdiagonalisation diagonalisation;
58 // Embedded objects
60 SortEigen<Field> _sort;
65 int mrhs;
67 // BLAS objects
69 int Nevec_acc; // Number of eigenvectors stored in the buffer evec_acc
70
71 void VectorPoly(std::vector<Field> &in,std::vector<Field> &out)
72 {
73 Field mrhs_in(mrhs_grid);
74 Field mrhs_out(mrhs_grid);
75 for(int r=0;r<in.size();r+=mrhs){
76 for(int rr=0;rr<mrhs;rr++){
77 int rrr = r+rr;
78 if(rrr >= in.size()) rrr = 0;
79 InsertSlice(in[rrr],mrhs_in,rr,0);
80 }
81 _poly(_Linop,mrhs_in,mrhs_out);
82 for(int rr=0;rr<mrhs;rr++){
83 int rrr = r+rr;
84 if(rrr < in.size()) {
85 ExtractSlice(out[rrr],mrhs_out,rr,0);
86 }
87 }
88 }
89 }
90 void SingleOperator(Field &in,Field &out)
91 {
92 Field mrhs_in(mrhs_grid);
93 Field mrhs_out(mrhs_grid);
94 for(int rr=0;rr<mrhs;rr++){
95 InsertSlice(in,mrhs_in,rr,0);
96 }
97 _Linop.HermOp(mrhs_in,mrhs_out);
98 ExtractSlice(out,mrhs_out,0,0);
99 }
100
101 // Constructor
103public:
105 GridBase * f_Grid,
106 GridBase * mrhs_Grid,
107 int _mrhs,
108 OperatorFunction<Field> & poly, // polynomial
109 int _Nstop, // really sought vecs
110 int _Nconv_test_interval, // conv check interval
111 int _Nu, // vecs in the unit block
112 int _Nk, // sought vecs
113 int _Nm, // total vecs
114 RealD _eresid, // resid in lmd deficit
115 int _MaxIter, // Max iterations
116 IRBLdiagonalisation _diagonalisation = IRBLdiagonaliseWithEigen)
117 : _Linop(Linop), _poly(poly),f_grid(f_Grid), mrhs_grid(mrhs_Grid),
118 Nstop(_Nstop), Nconv_test_interval(_Nconv_test_interval), mrhs(_mrhs),
119 Nu(_Nu), Nk(_Nk), Nm(_Nm),
120 Nblock_m(_Nm/_Nu), Nblock_k(_Nk/_Nu),
121 eresid(_eresid), MaxIter(_MaxIter),
122 diagonalisation(_diagonalisation),
123 Nevec_acc(_Nu)
124 { assert( (Nk%Nu==0) && (Nm%Nu==0) ); };
125
127 // Helpers
129 static RealD normalize(Field& v, int if_print=0)
130 {
131 RealD nn = norm2(v);
132 nn = sqrt(nn);
133 v = v * (1.0/nn);
134 return nn;
135 }
136
137 void orthogonalize(Field& w, std::vector<Field>& evec, int k, int if_print=0)
138 {
139 typedef typename Field::scalar_type MyComplex;
140 ComplexD ip;
141
142 for(int j=0; j<k; ++j){
143 ip = innerProduct(evec[j],w);
144 if(if_print)
145 if( norm(ip)/norm2(w) > 1e-14)
146 Glog<<"orthogonalize before: "<<j<<" of "<<k<<" "<< ip <<std::endl;
147 w = w - ip * evec[j];
148 if(if_print) {
149 ip = innerProduct(evec[j],w);
150 if( norm(ip)/norm2(w) > 1e-14)
151 Glog<<"orthogonalize after: "<<j<<" of "<<k<<" "<< ip <<std::endl;
152 }
153 }
154 assert(normalize(w,if_print) != 0);
155 }
156 void reorthogonalize(Field& w, std::vector<Field>& evec, int k)
157 {
158 orthogonalize(w, evec, k,1);
159 }
160
161 void orthogonalize(std::vector<Field>& w, int _Nu, std::vector<Field>& evec, int k, int if_print=0)
162 {
163 typedef typename Field::scalar_type MyComplex;
164 MyComplex ip;
165
166 for(int j=0; j<k; ++j){
167 for(int i=0; i<_Nu; ++i){
168 ip = innerProduct(evec[j],w[i]);
169 w[i] = w[i] - ip * evec[j];
170 }}
171 for(int i=0; i<_Nu; ++i)
172 assert(normalize(w[i],if_print) !=0);
173 }
174
175 void orthogonalize_blockhead(Field& w, std::vector<Field>& evec, int k, int Nu)
176 {
177 typedef typename Field::scalar_type MyComplex;
178 MyComplex ip;
179
180 for(int j=0; j<k; ++j){
181 ip = innerProduct(evec[j*Nu],w);
182 w = w - ip * evec[j*Nu];
183 }
184 normalize(w);
185 }
186
187 void calc(std::vector<RealD>& eval,
188 std::vector<Field>& evec,
189 const std::vector<Field>& src, int& Nconv, LanczosType Impl)
190 {
191 switch (Impl) {
192 case LanczosType::irbl:
193 calc_irbl(eval,evec,src,Nconv);
194 break;
195
196 case LanczosType::rbl:
197 calc_rbl(eval,evec,src,Nconv);
198 break;
199 }
200 }
201
202 void calc_irbl(std::vector<RealD>& eval,
203 std::vector<Field>& evec,
204 const std::vector<Field>& src, int& Nconv)
205 {
206 std::string fname = std::string(cname+"::calc_irbl()");
207 GridBase *grid = evec[0].Grid();
208 assert(grid == src[0].Grid());
209 assert( Nu = src.size() );
210
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;
219 if ( diagonalisation == IRBLdiagonaliseWithEigen ) {
220 Glog << "Diagonalisation is Eigen "<< std::endl;
221#ifdef USE_LAPACK
222 } else if ( diagonalisation == IRBLdiagonaliseWithLAPACK ) {
223 Glog << "Diagonalisation is LAPACK "<< std::endl;
224#endif
225 } else {
226 abort();
227 }
228 Glog << std::string(74,'*') << std::endl;
229
230 assert(Nm == evec.size() && Nm == eval.size());
231
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);
238
239 Eigen::MatrixXcd Qt = Eigen::MatrixXcd::Zero(Nm,Nm);
240 Eigen::MatrixXcd Q = Eigen::MatrixXcd::Zero(Nm,Nm);
241
242 std::vector<int> Iconv(Nm);
243 std::vector<Field> B(Nm,grid); // waste of space replicating
244
245 std::vector<Field> f(Nu,grid);
246 std::vector<Field> f_copy(Nu,grid);
247 Field v(grid);
248
249 Nconv = 0;
250
251 RealD beta_k;
252
253 // set initial vector
254 for (int i=0; i<Nu; ++i) {
255 Glog << "norm2(src[" << i << "])= "<< norm2(src[i]) << std::endl;
256 evec[i] = src[i];
257 orthogonalize(evec[i],evec,i);
258 // Glog << "norm2(evec[" << i << "])= "<< norm2(evec[i]) << std::endl;
259 }
260
261 // initial Nblock_k steps
262 for(int b=0; b<Nblock_k; ++b) blockwiseStep(lmd,lme,evec,f,f_copy,b);
263
264 // restarting loop begins
265 int iter;
266 for(iter = 0; iter<MaxIter; ++iter){
267
268 Glog <<"#Restart iteration = "<< iter << std::endl;
269 // additional (Nblock_m - Nblock_k) steps
270 for(int b=Nblock_k; b<Nblock_m; ++b) blockwiseStep(lmd,lme,evec,f,f_copy,b);
271
272 // getting eigenvalues
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];
277 }
278 }
279 Qt = Eigen::MatrixXcd::Identity(Nm,Nm);
280 diagonalize(eval2,lmd2,lme2,Nu,Nm,Nm,Qt,grid);
281 _sort.push(eval2,Nm);
282 // Glog << "#Ritz value before shift: "<< std::endl;
283 for(int i=0; i<Nm; ++i){
284 // std::cout.precision(13);
285 // std::cout << "[" << std::setw(4)<< std::setiosflags(std::ios_base::right) <<i<<"] ";
286 // std::cout << "Rval = "<<std::setw(20)<< std::setiosflags(std::ios_base::left)<< eval2[i] << std::endl;
287 }
288
289 //----------------------------------------------------------------------
290 if ( Nm>Nk ) {
291 // Glog <<" #Apply shifted QR transformations "<<std::endl;
292 //int k2 = Nk+Nu;
293 int k2 = Nk;
294
295 Eigen::MatrixXcd BTDM = Eigen::MatrixXcd::Identity(Nm,Nm);
296 Q = Eigen::MatrixXcd::Identity(Nm,Nm);
297
299
300 for(int ip=Nk; ip<Nm; ++ip){
301 Glog << " ip "<<ip<<" / "<<Nm<<std::endl;
302 shiftedQRDecompEigen(BTDM,Nu,Nm,eval2[ip],Q);
303 }
304
306
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);
312 }
313 }
314 for(int i=0; i<k2; ++i) evec[i] = B[i];
315
316 // reconstruct initial vector for additional pole space
317 blockwiseStep(lmd,lme,evec,f,f_copy,Nblock_k-1);
318
319 // getting eigenvalues
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];
324 }
325 }
326 Qt = Eigen::MatrixXcd::Identity(Nm,Nm);
327 diagonalize(eval2,lmd2,lme2,Nu,Nk,Nm,Qt,grid);
328 _sort.push(eval2,Nk);
329 // Glog << "#Ritz value after shift: "<< std::endl;
330 for(int i=0; i<Nk; ++i){
331 // std::cout.precision(13);
332 // std::cout << "[" << std::setw(4)<< std::setiosflags(std::ios_base::right) <<i<<"] ";
333 // std::cout << "Rval = "<<std::setw(20)<< std::setiosflags(std::ios_base::left)<< eval2[i] << std::endl;
334 }
335 }
336 //----------------------------------------------------------------------
337
338 // Convergence test
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);
345 }
346 }
347
348 Nconv = 0;
349 for(int i=0; i<Nk; ++i){
350
351 // _Linop.HermOp(B[i],v);
352 SingleOperator(B[i],v);
353
354 RealD vnum = real(innerProduct(B[i],v)); // HermOp.
355 RealD vden = norm2(B[i]);
356 eval2[i] = vnum/vden;
357 v -= eval2[i]*B[i];
358 RealD vv = norm2(v);
359 resid[i] = vv;
360
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;
365
366 // change the criteria as evals are supposed to be sorted, all evals smaller(larger) than Nstop should have converged
367 //if( (vv<eresid*eresid) && (i == Nconv) ){
368 if (vv<eresid*eresid) {
369 Iconv[Nconv] = i;
370 ++Nconv;
371 }
372
373 } // i-loop end
374
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;
381 }
382
383 if ( Nconv>=Nstop ) break;
384
385 } // end of iter loop
386
387 Glog << std::string(74,'*') << std::endl;
388 if ( Nconv<Nstop ) {
389 Glog << fname + " NOT converged ; Summary :\n";
390 } else {
391 Glog << fname + " CONVERGED ; Summary :\n";
392 // Sort convered eigenpairs.
393 eval.resize(Nconv);
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]];
398 }
399 _sort.push(eval,evec,Nconv);
400 }
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;
406
407 }
408
409
410 void calc_rbl(std::vector<RealD>& eval,
411 std::vector<Field>& evec,
412 const std::vector<Field>& src, int& Nconv)
413 {
414 std::string fname = std::string(cname+"::calc_rbl()");
415 GridBase *grid = evec[0].Grid();
416 assert(grid == src[0].Grid());
417 assert( Nu = src.size() );
418
419 int Np = (Nm-Nk);
420 if (Np > 0 && MaxIter > 1) Np /= MaxIter;
421 int Nblock_p = Np/Nu;
422 for(int i=0;i< evec.size();i++) evec[0].Advise()=AdviseInfrequentUse;
423
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;
433 if ( diagonalisation == IRBLdiagonaliseWithEigen ) {
434 Glog << "Diagonalisation is Eigen "<< std::endl;
435#ifdef USE_LAPACK
436 } else if ( diagonalisation == IRBLdiagonaliseWithLAPACK ) {
437 Glog << "Diagonalisation is LAPACK "<< std::endl;
438#endif
439 } else {
440 abort();
441 }
442 Glog << std::string(74,'*') << std::endl;
443
444 assert(Nm == evec.size() && Nm == eval.size());
445
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);
452
453 Eigen::MatrixXcd Qt = Eigen::MatrixXcd::Zero(Nm,Nm);
454 Eigen::MatrixXcd Q = Eigen::MatrixXcd::Zero(Nm,Nm);
455
456 std::vector<int> Iconv(Nm);
457// int Ntest=Nu;
458// std::vector<Field> B(Nm,grid); // waste of space replicating
459 std::vector<Field> B(1,grid); // waste of space replicating
460
461 std::vector<Field> f(Nu,grid);
462 std::vector<Field> f_copy(Nu,grid);
463 Field v(grid);
464
465 Nconv = 0;
466
467// RealD beta_k;
468
469 // set initial vector
470 for (int i=0; i<Nu; ++i) {
471 Glog << "norm2(src[" << i << "])= "<< norm2(src[i]) << std::endl;
472 evec[i] = src[i];
473 orthogonalize(evec[i],evec,i);
474 Glog << "norm2(evec[" << i << "])= "<< norm2(evec[i]) << std::endl;
475 }
476// exit(-43);
477
478 // initial Nblock_k steps
479 for(int b=0; b<Nblock_k; ++b) blockwiseStep(lmd,lme,evec,f,f_copy,b);
480
481 // restarting loop begins
482 int iter;
483 int Nblock_l, Nblock_r;
484 int Nl, Nr;
485 int Nconv_guess = 0;
486
487 for(iter = 0; iter<MaxIter; ++iter){
488
489 Glog <<"#Restart iteration = "<< iter << std::endl;
490
491 Nblock_l = Nblock_k + iter*Nblock_p;
492 Nblock_r = Nblock_l + Nblock_p;
493 Nl = Nblock_l*Nu;
494 Nr = Nblock_r*Nu;
495 eval2.resize(Nr);
496
497 // additional Nblock_p steps
498 for(int b=Nblock_l; b<Nblock_r; ++b) blockwiseStep(lmd,lme,evec,f,f_copy,b);
499
500 // getting eigenvalues
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];
505 }
506 }
507 Qt = Eigen::MatrixXcd::Identity(Nr,Nr);
508 diagonalize(eval2,lmd2,lme2,Nu,Nr,Nr,Qt,grid);
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;
515 }
516
517 // Convergence test
518 Glog <<" #Convergence test: "<<std::endl;
519 Nconv = 0;
520 for(int j = 0; j<Nr; j+=Nconv_test_interval){
521 B[0]=0.0;
522 if ( j/Nconv_test_interval == Nconv ) {
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);
529 }
530
531 SingleOperator(B[0],v);
532 // _Linop.HermOp(B[0],v);
533 RealD vnum = real(innerProduct(B[0],v)); // HermOp.
534 RealD vden = norm2(B[0]);
535 eval2[j] = vnum/vden;
536 v -= eval2[j]*B[0];
537 RealD vv = norm2(v);
538 resid[j] = vv;
539
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;
544
545 // change the criteria as evals are supposed to be sorted, all evals smaller(larger) than Nstop should have converged
546 //if( (vv<eresid*eresid) && (i == Nconv) ){
547 if (vv<eresid*eresid) {
548 Iconv[Nconv] = j;
549 ++Nconv;
550 }
551 } else {
552 break;
553 }
554 } // j-loop end
555
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;
562 }
563
564 (Nconv > 0 ) ? Nconv_guess = 1 + (Nconv-1)*Nconv_test_interval : Nconv_guess = 0;
565 if ( Nconv_guess >= Nstop ) break;
566
567 } // end of iter loop
568
569 Glog << std::string(74,'*') << std::endl;
570 if ( Nconv_guess < Nstop ) {
571 Glog << fname + " NOT converged ; Summary :\n";
572 } else {
573 Glog << fname + " CONVERGED ; Summary :\n";
574 Nstop = Nconv_guess; // Just take them all
575 // Sort convered eigenpairs.
576 std::vector<Field> Btmp(Nstop,grid); // waste of space replicating
577
578 for(int i=0; i<Nstop; ++i){
579 Btmp[i]=0.;
580 for(int k = 0; k<Nr; ++k){
581 Btmp[i].Checkerboard() = evec[k].Checkerboard();
582 Btmp[i] += evec[k]*Qt(k,i);
583 }
584 SingleOperator(Btmp[i],v);
585 // _Linop.HermOp(Btmp[i],v);
586 RealD vnum = real(innerProduct(Btmp[i],v)); // HermOp.
587 RealD vden = norm2(Btmp[i]);
588// eval2[j] = vnum/vden;
589 v -= vnum/vden*Btmp[i];
590 RealD vv = norm2(v);
591// resid[j] = vv;
592
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;
597 eval[i] = vnum/vden;
598 }
599 for(int i=0; i<Nstop; ++i) evec[i] = Btmp[i];
600 eval.resize(Nstop);
601 evec.resize(Nstop,grid);
602 _sort.push(eval,evec,Nstop);
603 }
604 Glog << std::string(74,'*') << std::endl;
605 Glog << " -- Iterations = "<< iter << "\n";
606 // Glog << " -- beta(k) = "<< beta_k << "\n";
607 Glog << " -- Nconv = "<< Nconv << "\n";
608 Glog << " -- Nconv (guess) = "<< Nconv_guess << "\n";
609 Glog << std::string(74,'*') << std::endl;
610
611 }
612
613private:
614 void blockwiseStep(std::vector<std::vector<ComplexD>>& lmd,
615 std::vector<std::vector<ComplexD>>& lme,
616 std::vector<Field>& evec,
617 std::vector<Field>& w,
618 std::vector<Field>& w_copy,
619 int b)
620 {
621 const RealD tiny = 1.0e-20;
622
623 int Nu = w.size();
624 int Nm = evec.size();
625 assert( b < Nm/Nu );
626
627 // converts block index to full indicies for an interval [L,R)
628 int L = Nu*b;
629 int R = Nu*(b+1);
630
631 Real beta;
632
633 assert((Nu%mrhs)==0);
634 std::vector<Field> in(mrhs,f_grid);
635 std::vector<Field> out(mrhs,f_grid);
636
637 // unnecessary copy. Can or should it be avoided?
638 int k_start = 0;
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];
642 VectorPoly(in,out);
643 for (int u=0; u<mrhs; ++u) w[k_start+u] = out[u];
644 // for (int u=0; u<mrhs; ++u) Glog << " out["<<u<<"] = "<<norm2(out[u])<<std::endl;
645 k_start +=mrhs;
646 }
647 // Glog << "LinAlg "<< std::endl;
648
649 if (b>0) {
650 for (int u=0; u<Nu; ++u) {
651 //for (int k=L-Nu; k<L; ++k) {
652 for (int k=L-Nu+u; k<L; ++k) {
653 w[u] = w[u] - evec[k] * conjugate(lme[u][k]);
654 }
655 }
656 }
657
658 // 4. αk:=(vk,wk)
659 //for (int u=0; u<Nu; ++u) {
660 // for (int k=L; k<R; ++k) {
661 // lmd[u][k] = innerProduct(evec[k],w[u]); // lmd = transpose of alpha
662 // }
663 // lmd[u][L+u] = real(lmd[u][L+u]); // force diagonal to be real
664 //}
665 for (int u=0; u<Nu; ++u) {
666 for (int k=L+u; k<R; ++k) {
667 lmd[u][k] = innerProduct(evec[k],w[u]); // lmd = transpose of alpha
668 // Glog <<"lmd "<<u<<" "<<k<<lmd[u][k] -conjugate(innerProduct(evec[u+L],w[k-L]))<<std::endl;
669 lmd[k-L][u+L] = conjugate(lmd[u][k]); // force hermicity
670 }
671 lmd[u][L+u] = real(lmd[u][L+u]); // force diagonal to be real
672 }
673
674 // 5. wk:=wk−αkvk
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];
678 }
679 w_copy[u] = w[u];
680 }
681 // Glog << "LinAlg done"<< std::endl;
682
683 // In block version, the steps 6 and 7 in Lanczos construction is
684 // replaced by the QR decomposition of new basis block.
685 // It results block version beta and orthonormal block basis.
686 // Here, QR decomposition is done by using Gram-Schmidt.
687 for (int u=0; u<Nu; ++u) {
688 for (int k=L; k<R; ++k) {
689 lme[u][k] = 0.0;
690 }
691 }
692
693 // re-orthogonalization for numerical stability
694 // Glog << "Gram Schmidt"<< std::endl;
695 orthogonalize(w,Nu,evec,R);
696 // QR part
697 for (int u=1; u<Nu; ++u) {
698 orthogonalize(w[u],w,u);
699 }
700 // Glog << "Gram Schmidt done "<< std::endl;
701
702 // Glog << "LinAlg "<< std::endl;
703 for (int u=0; u<Nu; ++u) {
704 //for (int v=0; v<Nu; ++v) {
705 for (int v=u; v<Nu; ++v) {
706 lme[u][L+v] = innerProduct(w[u],w_copy[v]);
707 }
708 lme[u][L+u] = real(lme[u][L+u]); // force diagonal to be real
709 }
710 //lme[0][L] = beta;
711
712 for (int u=0; u<Nu; ++u) {
713 // Glog << "norm2(w[" << u << "])= "<< norm2(w[u]) << std::endl;
714 assert (!isnan(norm2(w[u])));
715 for (int k=L+u; k<R; ++k) {
716 // Glog <<" In block "<< b << "," <<" beta[" << u << "," << k-L << "] = " << lme[u][k] << std::endl;
717 }
718 }
719 // Glog << "LinAlg done "<< std::endl;
720
721 if (b < Nm/Nu-1) {
722 for (int u=0; u<Nu; ++u) {
723 evec[R+u] = w[u];
724 }
725 }
726
727 }
728
729
730 void diagonalize_Eigen(std::vector<RealD>& eval,
731 std::vector<std::vector<ComplexD>>& lmd,
732 std::vector<std::vector<ComplexD>>& lme,
733 int Nu, int Nk, int Nm,
734 Eigen::MatrixXcd & Qt, // Nm x Nm
735 GridBase *grid)
736 {
737 assert( Nk%Nu == 0 && Nm%Nu == 0 );
738 assert( Nk <= Nm );
739 Eigen::MatrixXcd BlockTriDiag = Eigen::MatrixXcd::Zero(Nk,Nk);
740
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];
744 }
745 }
746
747 for ( int u=0; u<Nu; ++u ) {
748 for (int k=Nu; k<Nk; ++k ) {
749 BlockTriDiag(k-Nu,u+(k/Nu)*Nu) = conjugate(lme[u][k-Nu]);
750 BlockTriDiag(u+(k/Nu)*Nu,k-Nu) = lme[u][k-Nu];
751 }
752 }
753 //std::cout << BlockTriDiag << std::endl;
754
755 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXcd> eigensolver(BlockTriDiag);
756
757 for (int i = 0; i < Nk; i++) {
758 eval[Nk-1-i] = eigensolver.eigenvalues()(i);
759 }
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);
763 //Qt(Nk-1-i,j) = eigensolver.eigenvectors()(i,j);
764 //Qt(i,j) = eigensolver.eigenvectors()(i,j);
765 }
766 }
767 }
768
769#ifdef USE_LAPACK
770 void diagonalize_lapack(std::vector<RealD>& eval,
771 std::vector<std::vector<ComplexD>>& lmd,
772 std::vector<std::vector<ComplexD>>& lme,
773 int Nu, int Nk, int Nm,
774 Eigen::MatrixXcd & Qt, // Nm x Nm
775 GridBase *grid)
776 {
777 Glog << "diagonalize_lapack: Nu= "<<Nu<<" Nk= "<<Nk<<" Nm= "<<std::endl;
778 assert( Nk%Nu == 0 && Nm%Nu == 0 );
779 assert( Nk <= Nm );
780 Eigen::MatrixXcd BlockTriDiag = Eigen::MatrixXcd::Zero(Nk,Nk);
781
782 for ( int u=0; u<Nu; ++u ) {
783 for (int k=0; k<Nk; ++k ) {
784 // Glog << "lmd "<<u<<" "<<k<<" "<<lmd[u][k] -conjugate(lmd[u][k])<<std::endl;
785 BlockTriDiag(k,u+(k/Nu)*Nu) = lmd[u][k];
786 }
787 }
788
789 for ( int u=0; u<Nu; ++u ) {
790 for (int k=Nu; k<Nk; ++k ) {
791// Glog << "lme "<<u<<" "<<k<<" "<<lme[u][k] -conjugate(lme[u][k])<<std::endl;
792 BlockTriDiag(k-Nu,u+(k/Nu)*Nu) = conjugate(lme[u][k-Nu]);
793 BlockTriDiag(u+(k/Nu)*Nu,k-Nu) = lme[u][k-Nu];
794 }
795 }
796 //std::cout << BlockTriDiag << std::endl;
797//#ifdef USE_LAPACK
798 const int size = Nm;
799 MKL_INT NN = Nk;
800// double evals_tmp[NN];
801// double evec_tmp[NN][NN];
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();
811 }
812 }
813 MKL_INT evals_found;
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];
819// double work[lwork];
820 MKL_Complex16 *work = (MKL_Complex16 *) malloc(lwork*sizeof(MKL_Complex16));
821 MKL_INT isuppz[2*NN];
822 char jobz = 'V'; // calculate evals & evecs
823 char range = 'I'; // calculate all evals
824 // char range = 'A'; // calculate all evals
825 char uplo = 'U'; // refer to upper half of original matrix
826 char compz = 'I'; // Compute eigenvectors of tridiagonal matrix
827 int ifail[NN];
828 MKL_INT info;
829 int total = grid->_Nprocessors;
830 int node = grid->_processor;
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);
834 if (iu > NN) iu=NN;
835 Glog << "node "<<node<<"il "<<il<<"iu "<<iu<<std::endl;
836 double tol = 0.0;
837 if (1) {
838 memset(evals_tmp,0,sizeof(double)*NN);
839 if ( il <= NN){
840 zheevr(&jobz, &range, &uplo, &NN,
841 DD, &NN,
842 &vl, &vu, &il, &iu, // these four are ignored if second parameteris 'A'
843 &tol, // tolerance
844 &evals_found, evals_tmp, (MKL_Complex16*)evec_tmp, &NN,
845 isuppz,
846 work, &lwork,
847 rwork, &lrwork,
848 iwork, &liwork,
849 &info);
850// (double*)EE,
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];
856 if (il>1) {
857 evec_tmp[(i-(il-1))*NN+j].imag=0.;
858 evec_tmp[(i-(il-1))*NN+j].real=0.;
859 }
860 }
861 }
862 }
863 {
864 grid->GlobalSumVector(evals_tmp,NN);
865 grid->GlobalSumVector((double*)evec_tmp,2*NN*NN);
866 }
867 }
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++) {
872// Qt(j,Nk-1-i) = eigensolver.eigenvectors()(j,i);
873 Qt(j,Nk-1-i)=std::complex<double>
874 ( evec_tmp[i*Nk+j].real,
875 evec_tmp[i*Nk+j].imag);
876// ( evec_tmp[(Nk-1-j)*Nk+Nk-1-i].real,
877// evec_tmp[(Nk-1-j)*Nk+Nk-1-i].imag);
878
879 }
880 }
881
882if (1){
883 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXcd> eigensolver(BlockTriDiag);
884
885 for (int i = 0; i < Nk; i++) {
886 Glog << "eval = "<<i<<" " <<eval[Nk-1-i] <<" "<< eigensolver.eigenvalues()(i) <<std::endl;
887// eval[Nk-1-i] = eigensolver.eigenvalues()(i);
888 }
889 for (int i = 0; i < Nk; i++) {
890 for (int j = 0; j < Nk; j++) {
891// Qt(j,Nk-1-i) = eigensolver.eigenvectors()(j,i);
892// Glog<<"Qt "<<j<<" "<<Nk-1-i<<" = " <<Qt(j,Nk-1-i) <<" "<<eigensolver.eigenvectors()(j,i) <<std::endl;
893 MKL_Complex16 tmp = evec_tmp[i*Nk+j];
894// Glog<<"Qt "<<j<<" "<<Nk-1-i<<" = " <<evec_tmp[(Nk-1-j)*Nk+Nk-1-i].real<<" "<<
895//evec_tmp[(Nk-1-j)*Nk+Nk-1-i].imag <<" "<<eigensolver.eigenvectors()(j,i) <<std::endl;
896 if ( (i<5)&& (j<5))
897 Glog<<"Qt "<<j<<" "<<Nk-1-i<<" = " << norm(Qt(j,Nk-1-i))<<" "<<
898 norm(eigensolver.eigenvectors()(j,i)) <<std::endl;
899 }
900 }
901}
902// exit(-43);
903
904 free (evals_tmp);
905 free (evec_tmp);
906 free (DD);
907 free (work);
908 }
909#endif
910
911
912 void diagonalize(std::vector<RealD>& eval,
913 std::vector<std::vector<ComplexD>>& lmd,
914 std::vector<std::vector<ComplexD>>& lme,
915 int Nu, int Nk, int Nm,
916 Eigen::MatrixXcd & Qt,
917 GridBase *grid)
918 {
919 Qt = Eigen::MatrixXcd::Identity(Nm,Nm);
920 if ( diagonalisation == IRBLdiagonaliseWithEigen ) {
921 diagonalize_Eigen(eval,lmd,lme,Nu,Nk,Nm,Qt,grid);
922#ifdef USE_LAPACK
923 } else if ( diagonalisation == IRBLdiagonaliseWithLAPACK ) {
924 diagonalize_lapack(eval,lmd,lme,Nu,Nk,Nm,Qt,grid);
925#endif
926 } else {
927 assert(0);
928 }
929 }
930
931
933 std::vector<std::vector<ComplexD>>& lmd,
934 std::vector<std::vector<ComplexD>>& lme,
935 int Nu, int Nb, int Nk, int Nm,
936 Eigen::MatrixXcd& M)
937 {
938 // Glog << "unpackHermitBlockTriDiagMatToEigen() begin" << '\n';
939 assert( Nk%Nu == 0 && Nm%Nu == 0 );
940 assert( Nk <= Nm );
941 M = Eigen::MatrixXcd::Zero(Nk,Nk);
942
943 // rearrange
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];
947 }
948 }
949
950 for ( int u=0; u<Nu; ++u ) {
951 for (int k=Nu; k<Nk; ++k ) {
952 M(k-Nu,u+(k/Nu)*Nu) = conjugate(lme[u][k-Nu]);
953 M(u+(k/Nu)*Nu,k-Nu) = lme[u][k-Nu];
954 }
955 }
956 // Glog << "unpackHermitBlockTriDiagMatToEigen() end" << std::endl;
957 }
958
959
961 std::vector<std::vector<ComplexD>>& lmd,
962 std::vector<std::vector<ComplexD>>& lme,
963 int Nu, int Nb, int Nk, int Nm,
964 Eigen::MatrixXcd& M)
965 {
966 // Glog << "packHermitBlockTriDiagMatfromEigen() begin" << '\n';
967 assert( Nk%Nu == 0 && Nm%Nu == 0 );
968 assert( Nk <= Nm );
969
970 // rearrange
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);
974 }
975 }
976
977 for ( int u=0; u<Nu; ++u ) {
978 for (int k=Nu; k<Nk; ++k ) {
979 lme[u][k-Nu] = M(u+(k/Nu)*Nu,k-Nu);
980 }
981 }
982 // Glog << "packHermitBlockTriDiagMatfromEigen() end" <<std::endl;
983 }
984
985
986 // assume the input matrix M is a band matrix
987 void shiftedQRDecompEigen(Eigen::MatrixXcd& M, int Nu, int Nm,
988 RealD Dsh,
989 Eigen::MatrixXcd& Qprod)
990 {
991 // Glog << "shiftedQRDecompEigen() begin" << '\n';
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);
995
996 Mtmp = M;
997 for (int i=0; i<Nm; ++i ) {
998 Mtmp(i,i) = M(i,i) - Dsh;
999 }
1000
1001 Eigen::HouseholderQR<Eigen::MatrixXcd> QRD(Mtmp);
1002 Q = QRD.householderQ();
1003 R = QRD.matrixQR(); // upper triangular part is the R matrix.
1004 // lower triangular part used to represent series
1005 // of Q sequence.
1006
1007 // Glog << "shiftedQRDecompEigen() Housholder & QR" << '\n';
1008 // equivalent operation of Qprod *= Q
1009 //M = Eigen::MatrixXcd::Zero(Nm,Nm);
1010
1011 //for (int i=0; i<Nm; ++i) {
1012 // for (int j=0; j<Nm-2*(Nu+1); ++j) {
1013 // for (int k=0; k<2*(Nu+1)+j; ++k) {
1014 // M(i,j) += Qprod(i,k)*Q(k,j);
1015 // }
1016 // }
1017 //}
1018 //for (int i=0; i<Nm; ++i) {
1019 // for (int j=Nm-2*(Nu+1); j<Nm; ++j) {
1020 // for (int k=0; k<Nm; ++k) {
1021 // M(i,j) += Qprod(i,k)*Q(k,j);
1022 // }
1023 // }
1024 //}
1025
1026 Mtmp = Eigen::MatrixXcd::Zero(Nm,Nm);
1027
1028 // Glog << "shiftedQRDecompEigen() Mtmp create" << '\n';
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);
1033 }
1034 }
1035 }
1036 // Glog << "shiftedQRDecompEigen() Mtmp loop1" << '\n';
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);
1041 }
1042 }
1043 }
1044 // Glog << "shiftedQRDecompEigen() Mtmp loop2" << '\n';
1045
1046 //static int ntimes = 2;
1047 //for (int j=0; j<Nm-(ntimes*Nu); ++j) {
1048 // for (int i=ntimes*Nu+j; i<Nm; ++i) {
1049 // Mtmp(i,j) = 0.0;
1050 // }
1051 //}
1052 //ntimes++;
1053
1054 Qprod = Mtmp;
1055
1056 // equivalent operation of M = Q.adjoint()*(M*Q)
1057 Mtmp = Eigen::MatrixXcd::Zero(Nm,Nm);
1058
1059 for (int a=0, i=0, kmax=0; a<Nu+1; ++a) {
1060 for (int j=0; j<Nm-a; ++j) {
1061 i = j+a;
1062 kmax = (Nu+1)+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);
1066 }
1067 Mtmp(j,i) = conj(Mtmp(i,j));
1068 }
1069 }
1070 // Glog << "shiftedQRDecompEigen() Mtmp loop3" << '\n';
1071
1072 for (int i=0; i<Nm; ++i) {
1073 Mtmp(i,i) = real(Mtmp(i,i)) + Dsh;
1074 }
1075
1076 // Glog << "shiftedQRDecompEigen() Mtmp loop4" << '\n';
1077 M = Mtmp;
1078
1079 //M = Q.adjoint()*(M*Q);
1080 //for (int i=0; i<Nm; ++i) {
1081 // for (int j=0; j<Nm; ++j) {
1082 // if (i==j) M(i,i) = real(M(i,i));
1083 // if (j>i) M(i,j) = conj(M(j,i));
1084 // if (i-j > Nu || j-i > Nu) M(i,j) = 0.;
1085 // }
1086 //}
1087
1088 // Glog << "shiftedQRDecompEigen() end" <<std::endl;
1089 }
1090
1092 {
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);
1097
1098 A(0,0) = 12.0;
1099 A(0,1) = -51.0;
1100 A(0,2) = 4.0;
1101 A(1,0) = 6.0;
1102 A(1,1) = 167.0;
1103 A(1,2) = -68.0;
1104 A(2,0) = -4.0;
1105 A(2,1) = 24.0;
1106 A(2,2) = -41.0;
1107
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';
1112 }
1113 }
1114 Glog << std::endl;
1115
1116 Eigen::ColPivHouseholderQR<Eigen::MatrixXd> QRD(A);
1117
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';
1122 }
1123 }
1124 Glog << std::endl;
1125
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';
1131 }
1132 }
1133 Glog << std::endl;
1134
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';
1140 }
1141 }
1142 Glog << std::endl;
1143
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';
1149 }
1150 }
1151 Glog << std::endl;
1152
1153 Glog << "matrixR" << std::endl;
1154 R = QRD.matrixR();
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';
1158 }
1159 }
1160 Glog << std::endl;
1161
1162 Glog << "rank = " << QRD.rank() << std::endl;
1163 Glog << "threshold = " << QRD.threshold() << std::endl;
1164
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';
1170 }
1171 }
1172 Glog << std::endl;
1173
1174
1175 Glog << "QR decomposition without column pivoting" << std::endl;
1176
1177 A(0,0) = 12.0;
1178 A(0,1) = -51.0;
1179 A(0,2) = 4.0;
1180 A(1,0) = 6.0;
1181 A(1,1) = 167.0;
1182 A(1,2) = -68.0;
1183 A(2,0) = -4.0;
1184 A(2,1) = 24.0;
1185 A(2,2) = -41.0;
1186
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';
1191 }
1192 }
1193 Glog << std::endl;
1194
1195 Eigen::HouseholderQR<Eigen::MatrixXd> QRDplain(A);
1196
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';
1202 }
1203 }
1204 Glog << std::endl;
1205
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';
1210 }
1211 }
1212 Glog << std::endl;
1213 }
1214
1215};
1216
1218#undef Glog
1219#undef USE_LAPACK
1220
#define conj(a, b, i)
accelerator_inline Grid_simd< S, V > sqrt(const Grid_simd< S, V > &r)
B
accelerator_inline sobj eval(const uint64_t ss, const sobj &arg)
Definition Lattice_ET.h:103
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)
@ AdviseInfrequentUse
#define NAMESPACE_BEGIN(A)
Definition Namespace.h:35
#define NAMESPACE_END(A)
Definition Namespace.h:36
RealF Real
Definition Simd.h:65
std::complex< RealD > ComplexD
Definition Simd.h:79
double RealD
Definition Simd.h:61
void GlobalSumVector(RealF *, int N)
void orthogonalize(Field &w, std::vector< Field > &evec, int k, int if_print=0)
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 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 orthogonalize_blockhead(Field &w, std::vector< Field > &evec, int k, int Nu)
void calc(std::vector< RealD > &eval, std::vector< Field > &evec, const std::vector< Field > &src, int &Nconv, LanczosType Impl)
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)