-
Notifications
You must be signed in to change notification settings - Fork 0
/
lanczos_sparse.f90
executable file
·1743 lines (1397 loc) · 56.7 KB
/
lanczos_sparse.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!> initialize lanczos vector with a constant vector
subroutine InitLanczosVector(vector)
use prec
use sparse
implicit none
real(dp) :: norm
type (WTParVec), pointer :: vector
call WTParVecSetSeqValue(vector)
norm= WTParVecNorm(vector)
vector= 1d0/norm*vector
return
endsubroutine InitLanczosVector
!> Lanczos algorithm
!> NumLczVectors : number of Lanczos vectors
!> Mdim : The dimension of the Ham when it becomes a dense matrix
!> input: initial vector; Matrix stored in PARCSR format
!> output : Alpha(n), Beta(n)
subroutine lanczos_seqsparse_cpu_z(NumLczVectors, NumLczVectors_out, Mdim, nnz, icsr, jcsr, acsr, InitialVector, Alpha, Betan)
use prec
use wmpi
use sparse
use para, only : stdout, eps12, OmegaNum
implicit none
!* inout variables
!> number of lanczos vectors
integer, intent(in) :: NumLczVectors
integer, intent(inout) :: NumLczVectors_out
integer, intent(in) :: Mdim
integer, intent(in) :: nnz
integer, intent(in) :: icsr(Mdim+1)
integer, intent(in) :: jcsr(nnz)
complex(dp), intent(in) :: InitialVector(Mdim)
complex(dp), intent(in) :: acsr(nnz)
complex(dp), intent(out) :: Alpha (NumLczVectors)
complex(dp), intent(out) :: Betan (NumLczVectors)
!* local variables
!* loop index
integer :: i
real(dp) :: time_start, time_end, time_start0
!* Lanczos vector
complex(dp), allocatable :: vec_0(:), vec_1(:), vec_2(:)
!* temp variables
real(dp) :: error
logical :: converged
complex(dp) :: t_z1, t_z2, zdotc
real(dp), allocatable :: dos_old(:)
!if (cpuid==0) write(stdout, '(2x,a)')'>> Lanczos algorithm'
allocate(dos_old(OmegaNum))
allocate(vec_0(Mdim), vec_1(Mdim), vec_2(Mdim))
dos_old= 0d0
vec_0= 0d0; vec_1= 0d0; vec_2= 0d0
alpha= 0d0
betan= 0d0
call now(time_start0)
time_start= time_start0
!** run Lanczos procedure to get a set of lanczos vectors <<!
!* initialize the first lanczos vector
vec_0= InitialVector
!* normalize \phi_0
!* <\phi_0|\phi_0>
t_z1= zdotc(Mdim, vec_0, 1, InitialVector, 1)
if (dble(t_z1)<eps12) t_z1= eps12 !< avoiding computational error
!> the norm of the initial Lanczos vector should not be zero
betan(1)= dsqrt(dble(t_z1))
!* vec_0= vec_0/betan(1)
vec_0= 1d0/betan(1)* vec_0
!* perform |vec(1,:)>=Ham |vec(0,:)>
#if defined (INTELMKL)
call mkl_zcsrgemv('N', Mdim, acsr, icsr, jcsr, vec_0, vec_1)
#endif
!* a_0= <phi_0|H|phi_0>
alpha(1) = zdotc(Mdim, vec_0, 1, vec_1, 1)
!* |phi_1>'=H|phi_0>-a_0|phi_0>
!call WTSeqVecAxpy(-alpha(1), vec_0, vec_1)
call zaxpy(mdim, -alpha(1), vec_0, 1, vec_1, 1)
betan(2) = zdotc(Mdim, vec_1, 1, vec_1, 1)
if (dble(betan(2))<eps12) betan(2)= eps12 !< avoiding computational error
betan(2) = dsqrt(dble(betan(2)))
!* rescale vec_1
vec_1= 1d0/betan(2)*vec_1
!* time measurement
call now(time_end)
do i= 2, NumLczVectors
if (Mdim>500000.and.mod(i, 100)==0.and.cpuid.eq.0) &
write(stdout, '(2a, i9, " /", i10, a, f10.1, "s", a, f10.1, "s")') &
' In lanczos_seqsparse_cpu_z ', ' i/NumLczVectors', i, NumLczVectors, &
' time elapsed: ', time_end-time_start0, &
' time left: ', (NumLczVectors-i+1d0)*(time_end-time_start)
call now(time_start)
#if defined (INTELMKL)
call mkl_zcsrgemv('N', Mdim, acsr, icsr, jcsr, vec_1, vec_2)
#endif
!* calculate alpha, betan
t_z1= zdotc(Mdim, vec_1, 1, vec_2, 1)
alpha(i)= dble(t_z1)
!* vec_2= vec_2- betan(i)*vec_0- alpha(i)*vec_1
!call WTSeqVecAxpy(-alpha(i), vec_1, vec_2)
call zaxpy(mdim, -alpha(i), vec_1, 1, vec_2, 1)
!call WTSeqVecAxpy(-betan(i), vec_0, vec_2)
call zaxpy(mdim, -betan(i), vec_0, 1, vec_2, 1)
if (i.eq.NumLczVectors)cycle
t_z2= zdotc(Mdim, vec_2, 1, vec_2, 1)
if (dble(t_z2)<eps12) t_z2= eps12 !< avoiding computational error
betan(i+1)= dsqrt(dble(t_z2))
vec_0= vec_1
vec_2= 1/betan(i+1)* vec_2
vec_1= vec_2
!> check convergence
converged= .false.
if (mod(i, 10).eq.0)then
call lanczos_converged(alpha, betan, NumLczVectors, i, dos_old,&
error, .true., converged)
endif
if (converged)then
NumLczVectors_out= i
if (cpuid==0) then
write(stdout, *)' >> Lanczos procedure is converged at iteration: ', i
write(stdout, '(a, E16.5)')' >> Estimated error is : ', error
endif
exit
endif
call now(time_end)
enddo ! iteration
1000 continue
if (.not.converged)then
NumLczVectors_out= NumLczVectors
if (cpuid==0) then
write(stdout, '(a, i16)')' >> Lanczos procedure is not converged at iteration: ', NumLczVectors
write(stdout, '(a, E16.5)')' >> Estimated error is : ', error
endif
endif
!if (cpuid==0) then
! write(stdout, '(a)') ' ilcz alpha betan'
! do i = 1, NumLczVectors
! write(stdout, '(i5,40f16.8)')i,dble(alpha(i)), dble(betan(i))
! enddo
!endif
return
end subroutine lanczos_seqsparse_cpu_z
!> Lanczos algorithm
!> NumLczVectors : number of Lanczos vectors
!> Mdim : The dimension of the Ham when it becomes a dense matrix
!> input: initial vector; Matrix stored in PARCSR format
!> output : Alpha(n), Beta(n)
subroutine lanczos_parsparse_cpu_z(NumLczVectors, Mdim, Ham, InitialVector, Alpha, Betan)
use prec
use sparse
use wmpi
use para, only : stdout, eps12
implicit none
!* inout variables
!> number of lanczos vectors
integer, intent(in) :: NumLczVectors, Mdim
type (WTParCSR), intent(in) :: Ham
type (WTParVec), intent(in) :: InitialVector
complex(dp), intent(out) :: Alpha (NumLczVectors)
complex(dp), intent(out) :: Betan (NumLczVectors)
!* local variables
!* loop index
integer :: i
!* Lanczos vector
type (WTParVec), pointer :: vec_0
type (WTParVec), pointer :: vec_1
type (WTParVec), pointer :: vec_2
!* temp variables
complex(dp) :: t_z1
complex(dp) :: t_z2
if (cpuid==0) write(stdout, '(2x,a)')'>> Lanczos algorithm'
!print *, '%=============================================================%'
!print *, '| Lanczos on cpu with sparse matrix |'
!print *, '%=============================================================%'
vec_0=> Null()
vec_1=> Null()
vec_2=> Null()
allocate(vec_0, vec_1, vec_2)
#if defined (MPI)
call WTParVecCreate(mpi_cmw, Mdim, vec_0)
call WTParVecCreate(mpi_cmw, Mdim, vec_1)
call WTParVecCreate(mpi_cmw, Mdim, vec_2)
#endif
call WTParVecInitialize(vec_0)
call WTParVecInitialize(vec_1)
call WTParVecInitialize(vec_2)
alpha= 0d0
betan= 0d0
!** run Lanczos procedure to get a set of lanczos vectors <<!
!* initialize the first lanczos vector
vec_0= InitialVector
!* normalize \phi_0
!* <\phi_0|\phi_0>
t_z1 = vec_0*InitialVector
if (dble(t_z1)<eps12) t_z1= eps12 !< avoiding computational error
!> the norm of the initial Lanczos vector should not be zero
betan(1)= dsqrt(dble(t_z1))
!* vec_0= vec_0/betan(1)
vec_0= 1d0/betan(1)* vec_0
!* perform |vec(1,:)>=Ham |vec(0,:)>
call WTParCSRMatVec(Ham, vec_0, vec_1)
!* a_0= <phi_0|H|phi_0>
alpha(1)= vec_0*vec_1
!* |phi_1>'=H|phi_0>-a_0|phi_0>
call WTParVecAxpy(-alpha(1), vec_0, vec_1)
betan(2) = vec_1*vec_1
if (dble(betan(2))<eps12) betan(2)= eps12 !< avoiding computational error
betan(2) = dsqrt(dble(betan(2)))
!* rescale vec_1
vec_1= 1d0/betan(2)*vec_1
do i= 2, NumLczVectors
call WTParCSRMatVec(Ham, vec_1, vec_2)
!* calculate alpha, betan
t_z1= vec_1*vec_2
alpha(i)= dble(t_z1)
!* vec_2= vec_2- betan(i)*vec_0- alpha(i)*vec_1
call WTParVecAxpy(-alpha(i), vec_1, vec_2)
call WTParVecAxpy(-betan(i), vec_0, vec_2)
if (i.eq.NumLczVectors)cycle
t_z2= vec_2*vec_2
if (dble(t_z2)<eps12) t_z2= eps12 !< avoiding computational error
betan(i+1)= dsqrt(dble(t_z2))
vec_0= vec_1
vec_2= 1/betan(i+1)* vec_2
vec_1= vec_2
enddo
1000 continue
if (cpuid==0) then
write(stdout, '(a)') ' ilcz alpha betan'
do i = 1, NumLczVectors
write(stdout, '(i5,40f16.8)')i,dble(alpha(i)), dble(betan(i))
enddo
endif
call WTParVecDestroy(vec_0)
call WTParVecDestroy(vec_1)
call WTParVecDestroy(vec_2)
nullify(vec_0, vec_1, vec_2)
return
end subroutine lanczos_parsparse_cpu_z
subroutine LandauLevel_B_dos_Lanczos
!> we calculate the the spectrum with given magnetic field strength indicated by Magp and Magq,
!> a serials of kpoints kpoints(3, NK) and a serial of energies Omega(OmegaNum).
!> the output is dos_B_omega(Nk, OmegaNum)
!> Bx and By are global variables, you have to set it up in the current subroutine
use prec
use sparse
use wmpi
use mt19937_64
use para, only : Magq, Num_Wann, Bx, By, zi, pi, eta_arc, E_arc, &
OmegaNum, OmegaMin, OmegaMax, Magp, stdout, Magp_min, Magp_max, &
outfileindex, Single_KPOINT_3D_DIRECT,splen,Is_Sparse_Hr, eV2Hartree, &
MagneticSuperProjectedArea,ijmax,NumLCZVecs, NumRandomConfs, Add_Zeeman_Field
implicit none
!> magnetic field strength, this number should compatiable with the magnetic supercell
!> which means B*AreaOfMagneticSupercell=2*pi
!> The size of magnetic supercell is controled by Nq.
!> here Magp should be integer from 1 to Nq
!> energy interval
!> OmegaNum is defined in the module.f90 and read from the input.dat or wt.in
real, allocatable :: omega(:)
real(dp) :: time_start, time_end, time_start0
!> spectrum calculated
real(dp), allocatable :: dos_B_omega(:, :, :), dos_B_omega_mpi(:, :, :)
real(dp), allocatable :: n_int(:, :)
integer :: i, Mdim, Nq, ie, ib, ierr, it, ieta, iter
integer :: nnzmax, nnz, NumLczVectors, Nmag, NumLczVectors_out
real(dp) :: k3(3), thetaj, B0, energy, B0Tesla_quantumflux_magsupcell, eta
!> magnetic field strength
real(dp), allocatable :: mag_Tesla(:), flux(:)
!> left vector and right vector
complex(dp), allocatable :: InitialVector(:), InitialVectors(:, :)
!> Lanczos hamiltonian
complex(dp), allocatable :: Alpha(:), Betan(:)
integer, allocatable :: icsr(:), jcsr(:), iwk(:)
complex(dp), allocatable :: acsr(:)
complex(dp) :: zdotc, norm, theta
real(dp) :: continued_fraction
logical :: term
integer :: NumberofEta, ie_Earc
real(dp), allocatable :: eta_array(:), n_Earc(:)
real(dp) :: Bx_in_Tesla, By_in_Tesla, Bz_in_Tesla
integer(8) :: iseed
real(dp) :: rseed0
NumberofEta=9
allocate(eta_array(NumberofEta))
allocate(n_Earc(NumberofEta))
eta_array=(/0.1d0, 0.2d0, 0.4d0, 0.8d0, 1.0d0, 2d0, 4d0, 8d0, 10d0/)
eta_array= eta_array*Eta_Arc
Nq= Magq
Nmag= Magp_max- Magp_min +1
Mdim= Num_Wann*Nq
nnzmax= Num_wann*(2*ijmax+1+2)*Mdim
if(Is_Sparse_Hr) nnzmax=splen*Nq
if (cpuid==0) then
write(stdout, '(a,i8)')' Magnetic supercell is Nq= ', Nq
write(stdout, '(a,i18)')' Hamiltonian matrix dimension Mdim= ', Mdim
write(stdout, '(a,i18)')' nnzmax for lanczos is ', nnzmax
write(stdout, '(a)')' Maximal integer for 32-digit is 2,147,483,647'
write(stdout, '(a, f20.1, a)')' Memory for dos matrix is ', &
(Nmag)/1024d0/1024d0*9*OmegaNum*16d0*2d0, ' MB'
write(stdout, '(a, f20.1, a)')' Memory for sparse Hamiltonian is ', &
nnzmax/1024d0/1024d0*(4d0+4d0+16d0), ' MB'
write(stdout, '(a, f20.1, a)')' Memory for Lanczos vectors is ', &
Mdim/1024d0/1024d0*(4d0+NumRandomConfs)*16d0, ' MB'
write(stdout, '(a, f20.1, a)')' Memory for DOS storage is ', &
Nmag*OmegaNum*NumberofEta*2d0/1024d0/1024d0*8d0, ' MB'
endif
NumLczVectors = NumLCZVecs
NumLczVectors_out = NumLCZVecs
if (NumLczVectors>Mdim) then
NumLczVectors = Mdim-1
endif
allocate(omega(OmegaNum), mag_Tesla(Nmag), flux(Nmag))
allocate(n_int(0:Omeganum, NumberofEta))
allocate(dos_B_omega(Nmag, OmegaNum, NumberofEta), dos_B_omega_mpi(Nmag, OmegaNum, NumberofEta))
dos_B_omega= 0d0; dos_B_omega_mpi= 0d0
!> energy
do ie=1, OmegaNum
omega(ie)= OmegaMin+ (OmegaMax-OmegaMin)* (ie-1d0)/dble(OmegaNum-1)
enddo ! ie
!> deal with the magnetic field
!> first transform the Bx By into B*Cos\theta, B*Sin\theta
if (abs(By)<1e-8) then
theta= 0d0
elseif (By>0) then
theta = atan(Bx/By)
else
theta = atan(Bx/By)+ pi
endif
if (abs(theta)<1e-8 .and. Bx<0 ) theta= theta+ pi
!> when B0=2*pi/Nq, the magnetic flux in the magnetic supercell is 2*pi
!> the magnetic flux in the primitive cell is 2*pi/(Nq)
B0= 2d0*pi/dble(Nq)
B0= abs(B0)
!> transform it into Tesla
!> B=2*pi*\phi_0/S0, where \phi_0 is the quantum flux, S0 is the projected area
!> \phi_0 = h/2e h=6.62607004*1E-34, e= 1.6*1E-19
!> B0Tesla is the magnetic field when the flux of the magnetic-supercell is 2*pi
B0Tesla_quantumflux_magsupcell= 6.62607004d0*1E-34/2d0/1.6021766208d0*1E19/MagneticSuperProjectedArea*1E20
if (cpuid==0) then
write(stdout, '(a, 2E18.8)')' Magnetic field B in Tesla that makes the flux through '
write(stdout, '(a, 2E18.8)')' the magnetic supercell to be one quantum flux : ', B0Tesla_quantumflux_magsupcell
endif
!> The flux in the unit cell changes from 0 to 2*pi*Nmag/Nq
do ib= Magp_min, Magp_max
flux(ib-Magp_min+1)= B0* ib
mag_Tesla(ib-Magp_min+1)= B0Tesla_quantumflux_magsupcell* ib
enddo
k3=Single_KPOINT_3D_DIRECT
allocate(InitialVector(Mdim))
allocate(InitialVectors(Mdim, NumRandomConfs))
allocate(Alpha(NumLczVectors), Betan(NumLczVectors))
allocate(icsr(nnzmax), jcsr(nnzmax), acsr(nnzmax))
allocate(iwk (Mdim+1))
!* Get the Hamiltonian
time_start= 0d0
time_start0= 0d0
call now(time_start0)
time_start= time_start0
time_end = time_start0
!> Fix the initial vectors for each magnetic flux
do it=1, NumRandomConfs
call now(rseed0)
iseed=cpuid*NumRandomConfs*3132+ it+123431+int(rseed0)
call init_genrand64(iseed)
do i=1, Mdim
!call random_number(harvest= thetaj)
thetaj= genrand64_real2()
InitialVectors(i, it)= exp(zi*2d0*pi*thetaj)
enddo
enddo
do iter=1+cpuid, Nmag*NumRandomConfs, num_cpu
!do ib=1+cpuid, Nmag, num_cpu
! do it= 1, NumRandomConfs
ib= (iter-1)/NumRandomConfs+ 1
it= (iter-1 - (ib-1)*NumRandomConfs)+ 1
if (cpuid.eq.0) &
write(stdout, '(2a, i9, " /", i10, a, f10.1, "s", a, f10.1, "s")') &
'In LandauLevel_B_dos_Lanczos ', ' iter/Nmag*NumRandomConfs ', iter, Nmag*NumRandomConfs, &
' time elapsed: ', time_end-time_start0, &
' time left: ', ((Nmag*NumRandomConfs-iter-1d0)/dble(num_cpu))*(time_end-time_start)
call now(time_start)
!> here we set the magnetic field along the first vector of the SURFACE card
Bx= -flux(ib)
By= 0d0
!> TODO: should be revised in future, here we only put magnetic field along z direction
Bx_in_Tesla=0d0; By_in_Tesla= 0d0; Bz_in_Tesla= mag_Tesla(ib)
Alpha= 0d0; Betan= 0d0
icsr=0; jcsr=0;acsr=0d0;iwk=0d0
InitialVector(:)= InitialVectors(:, it)
norm= zdotc(Mdim, InitialVector, 1, InitialVector, 1)
InitialVector= InitialVector/dsqrt(dble(norm))
nnz= nnzmax
if(Is_Sparse_Hr) then
!> this subroutine will change the tight binding model according to the magnetic field strength
if (Add_Zeeman_Field) call add_zeeman_sparse_hr(Bx_in_Tesla, By_in_Tesla, Bz_in_Tesla)
!> get the Hamiltonian for the magnetic supercell whose size is Nq.
call ham_3Dlandau_sparseHR(nnz,Mdim,Nq,k3,acsr,jcsr,icsr)
else
call ham_3Dlandau_sparse1(nnz, Mdim, Nq, k3, acsr,jcsr,icsr)
end if
!> transform coo format to csr format
call ConvertCooToCsr(Mdim, nnz, acsr, icsr, jcsr, iwk)
call csr_sort_indices(Mdim, nnz, icsr, jcsr, acsr)
!> eleminate the same entries in the sparse matrix
call csr_sum_duplicates(Mdim, nnz, icsr, jcsr, acsr)
call csr_sort_indices(Mdim, nnz, icsr, jcsr, acsr)
call now(time_end)
if (cpuid.eq.0) write(stdout, '(a, f10.1, "s")') &
' Hamiltonian construction time cost :', time_end-time_start
!* doing lanczos procedure in order to get alpha, beta
call lanczos_seqsparse_cpu_z(NumLczVectors, NumLczVectors_out, Mdim, nnz, icsr, jcsr, acsr, InitialVector, Alpha, Betan)
!> Betan(1) is meaningless.
term = .True.
do ieta=1, NumberofEta
eta= eta_array(ieta)
do ie=1, OmegaNum
energy= omega(ie)
dos_B_omega(ib, ie, ieta)= dos_B_omega(ib, ie, ieta)+ &
continued_fraction(Alpha, Betan, energy, eta, NumLczVectors_out, term)
enddo
enddo
call now(time_end)
!enddo ! it= 1, NumRandomConfs
!enddo ! ib magnetic field
enddo ! combing ib and it
#if defined (MPI)
call mpi_allreduce(dos_B_omega, dos_B_omega_mpi, size(dos_B_omega), &
mpi_dp, mpi_sum, mpi_cmw, ierr)
#else
dos_B_omega_mpi= dos_B_omega
#endif
dos_B_omega_mpi= dos_B_omega_mpi/NumRandomConfs
outfileindex= outfileindex+ 1
if (cpuid.eq.0)then
open (unit=outfileindex, file='LandauLevel_B_dos.dat')
open (unit=outfileindex+1,file='wannierdiagram.dat')
write(outfileindex, '("#", a)')'Hofstadter butterfly (Landau level in (E,B) axis) '
write(outfileindex, '("#", a, i6)')'Magnetic supercell size : ', Magq
write(outfileindex+1, '("#", a)')'Wannier diagram '
write(outfileindex+1, '("#", a, i6)')'Magnetic supercell size : ', Magq
write(outfileindex, '("#", a, 3f12.6, a)')'k points in fractional coordinates (', k3, ')'
write(outfileindex, '("# Column", I5, 100I16)')(i, i=1, 12)
write(outfileindex, '("#", a15, 5a16)', advance='NO')'Flux', 'B(Tesla)', 'E(eV)'
write(outfileindex+1, '("# Column", I5, 100I16)')(i, i=1, 20)
write(outfileindex+1, '("#", a15, 2a16)', advance='NO')'Flux', 'B(Tesla)'
do ieta=1, NumberofEta-1
write(outfileindex, '(a16)', advance='NO') 'A(B, E)'
write(outfileindex+1, '(2a16)', advance='NO') 'n ', 'A(n, E)'
enddo
write(outfileindex, '(a16)') 'A(B, E)'
write(outfileindex+1, '(2a16)') 'n ', 'A(n, E)'
write(outfileindex, '("#", a, 20X, 300f16.2)')'Broadening \eta (meV): ', Eta_array(:)*1000d0/eV2Hartree
write(outfileindex+1, '("#", a, 5X, f26.2,300f32.2)')'Broadening \eta (meV): ', Eta_array(:)*1000d0/eV2Hartree
!> find n_int at EF
do ie=1, omeganum
if (omega(ie)>E_arc) then
ie_Earc= ie- 1
exit
endif
enddo
do ib=1, Nmag
n_int=0
do ie=1, omeganum
write(outfileindex, '(300f16.6)')flux(ib)/2d0/pi, mag_Tesla(ib), omega(ie)/eV2Hartree, dos_B_omega_mpi(ib, ie, :)
enddo
!> get number of electrons between the lowest energy level and E_arc
n_Earc= 0d0
do ie=1, ie_Earc
n_Earc(:)= n_Earc(:)+ dos_B_omega_mpi(ib, ie, :)
enddo
!> get number of electrons between the lowest energy level and omega(ie)
do ie=1, omeganum
n_int(ie, :)=n_int(ie-1, :)+ dos_B_omega_mpi(ib, ie, :)
enddo
!> set n(E)= n_int- n_Earc= \int_Earc^E \rho(\epsilon)d\epsilon
do ie=1, omeganum
n_int(ie, :)= n_int(ie, :)-n_Earc(:)
enddo
!do ieta=1, NumberofEta
! n_int(:, ieta)=n_int(:, ieta)/n_int(omeganum, ieta)
!enddo
do ie=1, omeganum
write(outfileindex+1,'(300f16.6)') flux(ib)/2d0/pi, mag_Tesla(ib), &
(n_int(ie, ieta), dos_B_omega_mpi(ib, ie, ieta), ieta=1, NumberofEta)
enddo
write(outfileindex+1, *) ' '
write(outfileindex, *) ' '
enddo
close(outfileindex)
write(stdout,*)'calculate Landau level spectrum in B-E mode successfully'
endif
outfileindex= outfileindex+ 2
if (cpuid.eq.0)then
open (unit=outfileindex, file='LandauLevel_B_dos.gnu')
write(outfileindex, '(a)')'#set terminal postscript enhanced color font ",24"'
write(outfileindex, '(a)')'set terminal pngcairo enhanced color font ",60" size 1920, 1680'
write(outfileindex, '(a)')"set output 'LandauLevel_B_dos.png'"
write(outfileindex, '(a)')'set pm3d'
write(outfileindex, '(a)')'set palette rgb 21,22,23'
write(outfileindex, '(a)')'#set isosamples 50,50'
write(outfileindex, '(a)')'set size 0.9,1'
write(outfileindex, '(a)')'set origin 0.05,0'
write(outfileindex, '(a)')'set view map'
write(outfileindex, '(a)')'unset ztics'
write(outfileindex, '(a)')'unset surface'
write(outfileindex, '(a)')'unset key'
write(outfileindex, '(a)')'#set xtics font ",24"'
write(outfileindex, '(a)')'#set ytics font ",24"'
write(outfileindex, '(a)')'#set xlabel font ",24"'
write(outfileindex, '(a)')'set ylabel "Energy (eV)"'
write(outfileindex, '(a, i6,a)')'set title "Hofstadter butterfly with Nq=', Nq, '" font ",40"'
write(outfileindex, '(a, f12.6, a, f12.6, a)')'set yrange [',OmegaMin/eV2Hartree,':',OmegaMax/eV2Hartree,']'
write(outfileindex, '(a)')'set xlabel "Phi per unit cell"'
write(outfileindex,*) '#set xlabel "{/Symbol F}/{/Symbol F}_0"'
write(outfileindex, '(a, f12.6, a)') 'set xrange [ 0.0000 :', maxval(flux/2d0/pi) ,']'
write(outfileindex, '(a)') "splot 'LandauLevel_B_dos.dat' u 1:3:(log($4)) w pm3d"
write(outfileindex, '(a)')'set xlabel "B (Tesla)"'
write(outfileindex, '(a, f12.6, a)') '#set xrange [ 0.0000 :', maxval(mag_Tesla) ,']'
write(outfileindex, '(a)') "#splot 'LandauLevel_B_dos.dat' u 2:3:(log($4)) w pm3d"
endif
outfileindex=outfileindex+1
if(cpuid == 0) then
open(unit=outfileindex, file='wannierdiagram.gnu')
write(outfileindex,*) 'set terminal pngcairo enhanced color font ",60" size 1920, 1680'
write(outfileindex,*) "set output 'wannierdiagram.png'"
write(outfileindex,*) 'set pm3d'
write(outfileindex, '(a)')'set size 0.9,1'
write(outfileindex, '(a)')'set origin 0.05,0'
write(outfileindex, '(a)')'set palette rgb 21,22,23'
write(outfileindex,*) '#set isosamples 50,50'
write(outfileindex,*) 'set view map'
write(outfileindex,*) 'unset ztics'
write(outfileindex,*) 'unset surface'
write(outfileindex,*) 'unset key'
write(outfileindex,*) 'set ylabel "n"'
write(outfileindex, '(a, i6, a)') 'set title "Wannier diagram with Nq=', Nq, '"'
write(outfileindex,*) '#set yrange [ ] noextend'
write(outfileindex,*) '#set xlabel "{/Symbol F}/{/Symbol F}_0"'
write(outfileindex,*) 'set xlabel "Phi/Phi_0 per unit cell"'
write(outfileindex,*) '#set xlabel "{/Symbol F}/{/Symbol F}_0"'
write(outfileindex,*) '#set xrange [ ] noextend'
write(outfileindex,*) "splot 'wannierdiagram.dat' u 1:3:(log($4)) w pm3d #lc palette"
end if
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
deallocate(InitialVector, Alpha, Betan, InitialVectors)
deallocate(icsr, jcsr, acsr, iwk)
return
end subroutine LandauLevel_B_dos_Lanczos
subroutine LandauLevel_k_dos_Lanczos
!> we calculate the the spectrum with given magnetic field strength indicated by Magp and Magq,
!> a serials of kpoints kpoints(3, NK) and a serial of energies Omega(OmegaNum).
!> the output is dos_k_omega(Nk, OmegaNum)
use prec
use sparse
use wmpi
use para, only : Magq, Num_Wann, Bx, By, zi, pi, eta_arc, Angstrom2atomic, &
OmegaNum, OmegaMin, OmegaMax, nk3_band, Magp, stdout, k3points, eV2Hartree, &
outfileindex, K3len_mag,splen,Is_Sparse_Hr,ijmax,NumLCZVecs, MagneticSuperProjectedArea, &
Nk3lines, k3line_mag_stop, k3line_name, NumRandomConfs
implicit none
!> magnetic field strength, this number should compatiable with the magnetic supercell
!> which means B*AreaOfMagneticSupercell=2*pi
!> The size of magnetic supercell is controled by Nq.
!> here Magp should be integer from 1 to Nq
!> energy interval
!> OmegaNum is defined in the module.f90 and read from the input.dat or wt.in
real, allocatable :: omega(:)
real(dp) :: time_start, time_end, time_start0
!> spectrum calculated
real(dp), allocatable :: dos_k_omega(:, :), dos_k_omega_mpi(:, :)
integer :: i, Mdim, Nq, ie, ik, ierr, it
integer :: nnzmax, nnz, NumLczVectors, NumLczVectors_out
real(dp) :: k3(3), thetaj, B0, energy, B0Tesla, B0Tesla_quantumflux_magsupcell
!> left vector and right vector
complex(dp), allocatable :: InitialVector(:)
!> Lanczos hamiltonian
complex(dp), allocatable :: Alpha(:), Betan(:)
integer, allocatable :: icsr(:), jcsr(:), iwk(:)
complex(dp), allocatable :: acsr(:)
complex(dp) :: zdotc, norm, theta
real(dp) :: continued_fraction
logical :: term
Nq= Magq
Mdim= Num_Wann*Magq
nnzmax= Num_wann*(2*ijmax+1)*Mdim
if(Is_Sparse_Hr) nnzmax=splen*Nq
!> need to be checked
!if(Is_Sparse_Hr) nnzmax=Nq*splen
NumLczVectors = NumLCZVecs
NumLczVectors_out = NumLCZVecs
if (NumLczVectors>Mdim) then
NumLczVectors = Mdim-1
endif
allocate(omega(OmegaNum))
allocate(dos_k_omega(nk3_band, OmegaNum), dos_k_omega_mpi(nk3_band, OmegaNum))
dos_k_omega= 0d0; dos_k_omega_mpi= 0d0
!> energy
do ie=1, OmegaNum
omega(ie)= OmegaMin+ (OmegaMax-OmegaMin)* (ie-1d0)/dble(OmegaNum-1)
enddo ! ie
!> deal with the magnetic field
!> first transform the Bx By into B*Cos\theta, B*Sin\theta
if (abs(By)<1e-8) then
theta= 0d0
elseif (By>0) then
theta = atan(Bx/By)
else
theta = atan(Bx/By)+ pi
endif
if (abs(theta)<1e-8 .and. Bx<0 ) theta= theta+ pi
!> the magnetic flux in the magnetic supercell is 2*pi*Magp
B0= 2d0*pi/dble(Nq)* Magp
B0= abs(B0)
Bx= B0* Cos(theta)
By= B0* Sin(theta)
!> transform it into Tesla
!> B=2*pi*\phi_0/S0, where \phi_0 is the quantum flux, S0 is the projected area
!> \phi_0 = h/2e h=6.62607004*1E-34, e= 1.6*1E-19
!> B0Tesla_quantumflux_magsupcell is the magnetic field that makes the flux through the magnetic
!> supercell to be one quantum flux
B0Tesla_quantumflux_magsupcell= 6.62607004*1E-34/2d0/1.6021766208*1E19/MagneticSuperProjectedArea*1E20
B0Tesla= 6.62607004*1E-34/2d0/1.6021766208*1E19/MagneticSuperProjectedArea*1E20*Magp
if (cpuid==0) then
write(stdout, '(a, 2E18.8)')' Magnetic field B in Tesla that makes the flux through '
write(stdout, '(a, 2E18.8)')' the magnetic supercell to be one quantum flux : ', B0Tesla_quantumflux_magsupcell
write(stdout, '(a, 2E18.8)')' Magnetic field B in Tesla that fits to Magp : ', B0Tesla
write(stdout, '(a, 2f18.8)')' Magnetic field Bx, By= ', Bx, By
endif
allocate(InitialVector(Mdim))
allocate(Alpha(NumLczVectors), Betan(NumLczVectors))
allocate(icsr(nnzmax), jcsr(nnzmax), acsr(nnzmax))
allocate(iwk (Mdim+1))
!* Get the Hamiltonian
time_start= 0d0
time_start0= 0d0
call now(time_start0)
time_start= time_start0
time_end = time_start0
do ik=1+cpuid, nk3_band, num_cpu
do it=1, NumRandomConfs
if (cpuid.eq.0) &
write(stdout, '(2a, i9, " /", i10, a, f10.1, "s", a, f10.1, "s")') &
'In LandauLevel_k_dos_Lanczos ', ' ik/NK ', ik, nk3_band, &
' time elapsed: ', time_end-time_start0, &
' time left: ', ((nk3_band-ik)/dble(num_cpu))*(time_end-time_start)*NumRandomConfs
call now(time_start)
Alpha= 0d0; Betan= 0d0
icsr=0; jcsr=0;acsr=0d0;iwk=0d0
do i=1, Mdim
call random_number(harvest= thetaj)
InitialVector(i)= exp(zi*2d0*pi*thetaj)
!InitialVector(i)= 1d0
enddo
norm= zdotc(Mdim, InitialVector, 1, InitialVector, 1)
InitialVector= InitialVector/dsqrt(dble(norm))
if (cpuid==0) write(stdout, '(a, 2i10)') 'LandauLevel_k_DOS', ik, nk3_band
k3= k3points(:, ik)
nnz= nnzmax
if(Is_Sparse_Hr) then
call ham_3Dlandau_sparseHR(nnz,Mdim,NQ,k3,acsr,jcsr,icsr)
else
call ham_3Dlandau_sparse1(nnz, Mdim, Nq, k3, acsr,jcsr,icsr)
end if
!> transform coo format to csr format
call ConvertCooToCsr(Mdim, nnz, acsr, icsr, jcsr, iwk)
call csr_sort_indices(Mdim, nnz, icsr, jcsr, acsr)
!> eleminate the same entries in the sparse matrix
call csr_sum_duplicates(Mdim, nnz, icsr, jcsr, acsr)
call csr_sort_indices(Mdim, nnz, icsr, jcsr, acsr)
!* doing lanczos procedure in order to get alpha, beta
call lanczos_seqsparse_cpu_z(NumLczVectors, NumLczVectors_out, Mdim, nnz, icsr, jcsr, acsr, InitialVector, Alpha, Betan)
!> Betan(1) is meaningless.
term = .False.
do ie=1, OmegaNum
energy= omega(ie)
dos_k_omega(ik, ie)=dos_k_omega(ik, ie)+ &
continued_fraction(Alpha,Betan,energy,eta_arc,NumLczVectors_out, term)
enddo
call now(time_end)
enddo ! it= 1, NumRandomConfs
enddo ! ik
#if defined (MPI)
call mpi_allreduce(dos_k_omega, dos_k_omega_mpi, size(dos_k_omega), &
mpi_dp, mpi_sum, mpi_cmw, ierr)
#else
dos_k_omega_mpi= dos_k_omega
#endif
dos_k_omega_mpi= dos_k_omega_mpi/NumRandomConfs
outfileindex= outfileindex+ 1
if (cpuid.eq.0)then
open (unit=outfileindex, file='LandauLevel_k_dos.dat')
write(outfileindex, '("#", a14, f16.8)')'Magnetic field in Tesla: ', B0Tesla
write(outfileindex, '("#", a14, 2a15, a)')'k', 'E(eV)', 'LDOS'
do ik=1, nk3_band
do ie=1, omeganum
write(outfileindex, '(30f16.8)')K3len_mag(ik)*Angstrom2atomic, omega(ie)/eV2Hartree, dos_k_omega_mpi(ik, ie)
enddo
write(outfileindex, *) ' '
enddo
close(outfileindex)
write(stdout,*)'calculate Landau level spectrum in k-E mode successfully'
endif
outfileindex= outfileindex+ 1
if (cpuid.eq.0)then
open (unit=outfileindex, file='LandauLevel_k_dos.gnu')
write(outfileindex, '(a)')'#set terminal postscript enhanced color font ",24"'
write(outfileindex, '(a)')'set terminal pngcairo enhanced color font ",70" size 1920, 1680'
write(outfileindex, '(a)')"set output 'LandauLevel_k_dos.png'"
write(outfileindex, '(a)')'set pm3d'
write(outfileindex, '(a)')'#set isosamples 50,50'
write(outfileindex, '(a)')'set size 0.85, 1'
write(outfileindex, '(a)')'set origin 0.1, 0'
write(outfileindex, '(a)')'set view map'
write(outfileindex, '(a)')'unset ztics'
write(outfileindex, '(a)')'unset surface'
write(outfileindex, '(a)')'unset key'
write(outfileindex, '(a)')'#set xtics font ",24"'
write(outfileindex, '(a)')'#set ytics font ",24"'
write(outfileindex, '(a)')'#set xlabel font ",24"'
write(outfileindex, '(a)')'set ylabel "Energy (eV)"'
write(outfileindex, '(a, i6,a)')'set title "Landau level with Nq=', Nq, '"'
write(outfileindex, '(a, f18.5, a)')'set xrange [0: ', maxval(K3len_mag*Angstrom2atomic), ']'
write(outfileindex, '(a, f18.5, a, f18.5, a)')'set yrange [', OmegaMin/eV2Hartree, ':', OmegaMax/eV2Hartree, ']'
write(outfileindex, 202, advance="no") (k3line_name(i), k3line_mag_stop(i), i=1, Nk3lines)
write(outfileindex, 203)k3line_name(nk3lines+1), k3line_mag_stop(Nk3lines+1)
do i=1, nk3lines-1
write(outfileindex, 204)k3line_mag_stop(i+1), OmegaMin, k3line_mag_stop(i+1), OmegaMax
enddo
write(outfileindex, '(a)')'set pm3d interpolate 2,2'
write(outfileindex, '(2a)')"splot 'LandauLevel_k_dos.dat' u 1:2:(log($3)) w pm3d"
close(outfileindex)
endif
202 format('set xtics (',:20('"',A1,'" ',F8.5,','))
203 format(A1,'" ',F8.5,')')
204 format('set arrow from ',F8.5,',',F10.5,' to ',F8.5,',',F10.5, ' nohead front lw 3')
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
deallocate(InitialVector, Alpha, Betan)
deallocate(icsr, jcsr, acsr, iwk)
return
end subroutine LandauLevel_k_dos_Lanczos
subroutine bulkbandk_dos_lanczos
!> we calculate the the spectrum with given magnetic field strength indicated by Magp and Magq,
!> a serials of kpoints kpoints(3, NK) and a serial of energies Omega(OmegaNum).
!> the output is dos_k_omega(Nk, OmegaNum)
use prec
use sparse
use wmpi
use para, only : Magq, Num_Wann, Bx, By, zi, pi, eta_arc, Angstrom2atomic, &
OmegaNum, OmegaMin, OmegaMax, nk3_band, Magp, stdout, k3points, &
outfileindex, K3len,splen,Is_Sparse_Hr,ijmax,NumLCZVecs, eV2Hartree
implicit none
!> magnetic field strength, this number should compatiable with the magnetic supercell
!> which means B*AreaOfMagneticSupercell=2*pi
!> The size of magnetic supercell is controled by Nq.
!> here Magp should be integer from 1 to Nq
!> energy interval
!> OmegaNum is defined in the module.f90 and read from the input.dat or wt.in
real, allocatable :: omega(:)
!> spectrum calculated
real(dp), allocatable :: dos_k_omega(:, :), dos_k_omega_mpi(:, :)
integer :: i, Mdim, Nq, ie, ik, ierr
integer :: nnzmax, nnz, NumLczVectors, NumLczVectors_out
real(dp) :: k3(3), thetaj, B0, energy
!> left vector and right vector
complex(dp), allocatable :: InitialVector(:)
!> Lanczos hamiltonian
complex(dp), allocatable :: Alpha(:), Betan(:)
integer, allocatable :: icsr(:), jcsr(:), iwk(:)
complex(dp), allocatable :: acsr(:)
complex(dp) :: zdotc, norm, theta
real(dp) :: continued_fraction
logical :: term
Nq= 1
Mdim= Num_Wann*Nq
nnzmax= Num_wann*(2*ijmax+1)*Mdim
if(Is_Sparse_Hr) nnzmax=splen*Nq
!> need to be checked
!if(Is_Sparse_Hr) nnzmax=Nq*splen
NumLczVectors = NumLCZVecs
if (NumLczVectors>Mdim) then
NumLczVectors = Mdim-1
endif
allocate(omega(OmegaNum))
allocate(dos_k_omega(nk3_band, OmegaNum), dos_k_omega_mpi(nk3_band, OmegaNum))
dos_k_omega= 0d0; dos_k_omega_mpi= 0d0
!> energy
do ie=1, OmegaNum
omega(ie)= OmegaMin+ (OmegaMax-OmegaMin)* (ie-1d0)/dble(OmegaNum-1)
enddo ! ie
!> deal with the magnetic field
!> first transform the Bx By into B*Cos\theta, B*Sin\theta
if (abs(By)<1e-8) then
theta= 0d0
elseif (By>0) then
theta = atan(Bx/By)
else
theta = atan(Bx/By)+ pi
endif
if (abs(theta)<1e-8 .and. Bx<0 ) theta= theta+ pi
!> the magnetic flux in the magnetic supercell is 2*pi*iq
allocate(InitialVector(Mdim))
allocate(Alpha(NumLczVectors), Betan(NumLczVectors))