-
Notifications
You must be signed in to change notification settings - Fork 0
/
landau_level_sparse.f90
executable file
·1788 lines (1471 loc) · 57.5 KB
/
landau_level_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
!include "mkl.fi"
!~ include "mkl_dss.f90"
subroutine ham_3Dlandau_sparse1(nnz, Ndimq, Nq, k, acoo,jcoo,icoo)
!> We generate the TB hamiltonian under magnetic field with dense formated HmnR.
!> At the begining, we generate coo-formated hamlitonian. However, in the end,
!> the coo-formated hamlitonian will be transform into csr format
use para
implicit none
integer, intent(inout) :: nnz
integer, intent(in) :: Ndimq
integer, intent(in) :: Nq
real(dp), intent(in) :: k(3)
! complex(dp) :: ham_landau(Ndimq, Ndimq)
complex(dp), intent(inout) :: acoo(nnz)
integer, intent(inout) :: jcoo(nnz)
integer, intent(inout) :: icoo(nnz)
!> inta-hopping for the supercell
complex(dp), allocatable :: H00(:, :)
!> inter-hopping for the supercell
complex(dp), allocatable :: H01(:, :)
complex(dp), allocatable :: H02(:, :)
! loop index
integer :: i1, i2,i, i_t, j_t
integer :: coo1,coo2,ncoo
! loop index
integer :: iR
! index used to sign irvec
real(dp) :: ia,ib,ic
integer :: ia1, ia2
integer :: istart1, istart2
integer :: iend1, iend2
integer :: inew_ic
!> nwann= Num_wann/2
integer :: nwann
integer, allocatable :: orbital_start(:)
!> new index used to sign irvec
real(dp) :: new_ia,new_ib,new_ic
!> wave vector k times lattice vector R
real(Dp) :: kdotr
real(dp) :: phase
complex(dp) :: ratio
complex(dp) :: fac,tmp
real(dp) :: Rp1(3)
real(dp) :: Rp2(3)
real(dp) :: R1(3)
real(dp) :: R2(3)
real(dp) :: Ri(3)
real(dp) :: Rj(3)
real(dp) :: tau1(3)
real(dp) :: tau2(3)
!> calculate phase induced by magnetic field
real(dp),external :: phase1,phase2
!> for leaking purpose, only will be allocated when nnz is not big enough
!complex(dp), allocatable :: acoo_extend(:)
!integer, allocatable :: icoo_extend(:)
!integer, allocatable :: jcoo_extend(:)
!> a work array for coocsr
integer(li), allocatable, target :: iwk (:)
allocate( H00( Num_wann, Num_wann))
allocate( H01( Num_wann, Num_wann))
allocate( H02( Num_wann, Num_wann))
H00= zzero
H01= zzero
H02= zzero
nwann= Num_wann/2
allocate( orbital_start(Origin_cell%Num_atoms+ 1))
orbital_start= 0
orbital_start(1)= 1
do ia1=1, Origin_cell%Num_atoms
orbital_start(ia1+1)= orbital_start(ia1)+ Origin_cell%nprojs(ia1)
enddo
! write(*,*) Origin_cell%Num_atoms
!> calculate intra-hopping
acoo=zzero
ncoo=0
! i1 column index, sweep magnetic supercells
do i1=1, Nq
! i2 row index, sweep magnetic supercells
do i2=1, Nq
H00=zzero
if (abs(i2-i1)> ijmax) cycle
!> sum over R points to get H(k1, k2)
do iR=1, Nrpts ! iR'th in lattice vector
ia=dble(irvec(1,iR))
ib=dble(irvec(2,iR))
ic=dble(irvec(3,iR))
!> new lattice
call latticetransform(ia, ib, ic, new_ia, new_ib, new_ic)
inew_ic= int(new_ic)
if (inew_ic /= (i2-i1)) cycle
!> exp(i k.R)
kdotr= k(1)*new_ia+ k(2)*new_ib
ratio= cos(2d0*pi*kdotr)+ zi*sin(2d0*pi*kdotr)
Rp1= (i1-1)*Ruc_new
Rp2= new_ia*Rua_new+ new_ib*Rub_new+ (i2-1)*Ruc_new
do ia1=1, Origin_cell%Num_atoms
do ia2=1, Origin_cell%Num_atoms
R1= Origin_cell%Atom_position_cart(:, ia1)
R2= Origin_cell%Atom_position_cart(:, ia2)
call rotate(R1, tau1)
call rotate(R2, tau2)
Ri= Rp1+ tau1
Rj= Rp2+ tau2
phase=phase2(Ri,Rj)
fac= cos(phase)+ zi*sin(phase)
istart1= orbital_start(ia1)
iend1= orbital_start(ia1+1)- 1
istart2= orbital_start(ia2)
iend2= orbital_start(ia2+1)- 1
H00 ( istart1:iend1, istart2:iend2) &
= H00 ( istart1:iend1, istart2:iend2) &
+ HmnR( istart1:iend1, istart2:iend2, iR)*ratio/ndegen(iR)* fac
!> there is soc term in the hr file
if (soc>0) then
istart1= orbital_start(ia1) + Nwann
iend1= orbital_start(ia1+1)- 1 + Nwann
istart2= orbital_start(ia2)
iend2= orbital_start(ia2+1)- 1
H00 ( istart1:iend1, istart2:iend2) &
= H00 ( istart1:iend1, istart2:iend2) &
+ HmnR( istart1:iend1, istart2:iend2, iR)*ratio/ndegen(iR)* fac
istart1= orbital_start(ia1)
iend1= orbital_start(ia1+1)- 1
istart2= orbital_start(ia2) + Nwann
iend2= orbital_start(ia2+1)- 1 + Nwann
H00 ( istart1:iend1, istart2:iend2) &
= H00 ( istart1:iend1, istart2:iend2) &
+ HmnR( istart1:iend1, istart2:iend2, iR)*ratio/ndegen(iR)* fac
istart1= orbital_start(ia1) + Nwann
iend1= orbital_start(ia1+1)- 1 + Nwann
istart2= orbital_start(ia2) + Nwann
iend2= orbital_start(ia2+1)- 1 + Nwann
H00 ( istart1:iend1, istart2:iend2) &
= H00 ( istart1:iend1, istart2:iend2) &
+ HmnR( istart1:iend1, istart2:iend2, iR)*ratio/ndegen(iR)* fac
endif ! soc
enddo ! ia2
enddo ! ia1
enddo ! iR
!> get the sparse matrix elements
do coo2=1,Num_wann
do coo1=1,Num_wann
tmp=H00(coo1,coo2)
if(abs(tmp)/eV2Hartree >eps6) then
j_t=coo1+(i1-1)*Num_wann
i_t=coo2+(i2-1)*Num_wann
!if (i_t>=j_t) then
ncoo=ncoo+1
jcoo(ncoo)=j_t
icoo(ncoo)=i_t
acoo(ncoo)=tmp !0.01
!endif
end if
end do
end do
H00=0d0
enddo ! i2
enddo ! i1
if (cpuid.eq.0) write(stdout,*) 'nnz, nnzmax after H00 : ', ncoo, nnz
if (nnz< ncoo) then
write(*, *)'Please increase nnz in the sparse.f90'
write(*,*) 'nnz, nnz after H00 : ', ncoo, nnz
stop
endif
!>> calculate inter-hopping between layers
! i1 column index
do i1=1, Nq
! i2 row index
do i2=1, Nq
H01=zzero
if (abs(i2+Nq -i1)> ijmax) cycle
!> sum over R points to get H(k1, k2)
do iR=1, Nrpts
ia=dble(irvec(1,iR))
ib=dble(irvec(2,iR))
ic=dble(irvec(3,iR))
!> new lattice
call latticetransform(ia, ib, ic, new_ia, new_ib, new_ic)
inew_ic= int(new_ic)
if (inew_ic /= (i2+ Nq -i1)) cycle
!> exp(i k.R)
kdotr= k(1)*new_ia+ k(2)*new_ib
ratio= cos(2d0*pi*kdotr)+ zi*sin(2d0*pi*kdotr)
Rp1= (i1-1)*Ruc_new
Rp2= new_ia*Rua_new+ new_ib*Rub_new+ (i2-1+ Nq)*Ruc_new
do ia1=1, Origin_cell%Num_atoms
do ia2=1, Origin_cell%Num_atoms
R1= Origin_cell%Atom_position_cart(:, ia1)
R2= Origin_cell%Atom_position_cart(:, ia2)
call rotate(R1, tau1)
call rotate(R2, tau2)
Ri= Rp1+ tau1
Rj= Rp2+ tau2
phase=phase2(Ri,Rj)
fac= cos(phase)+ zi*sin(phase)
istart1=orbital_start(ia1)
iend1= orbital_start(ia1+1)- 1
istart2= orbital_start(ia2)
iend2= orbital_start(ia2+1)- 1
H01 ( istart1:iend1, istart2:iend2) &
= H01 ( istart1:iend1, istart2:iend2) &
+ HmnR( istart1:iend1, istart2:iend2, iR)*ratio/ndegen(iR)* fac
!> there is soc term in the hr file
if (soc>0) then
istart1= orbital_start(ia1) + Nwann
iend1= orbital_start(ia1+1)- 1 + Nwann
istart2= orbital_start(ia2)
iend2= orbital_start(ia2+1)- 1
H01 ( istart1:iend1, istart2:iend2) &
= H01 ( istart1:iend1, istart2:iend2) &
+ HmnR( istart1:iend1, istart2:iend2, iR)*ratio/ndegen(iR)* fac
istart1= orbital_start(ia1)
iend1= orbital_start(ia1+1)- 1
istart2= orbital_start(ia2) + Nwann
iend2= orbital_start(ia2+1)- 1 + Nwann
H01 ( istart1:iend1, istart2:iend2) &
= H01 ( istart1:iend1, istart2:iend2) &
+ HmnR( istart1:iend1, istart2:iend2, iR)*ratio/ndegen(iR)* fac
istart1= orbital_start(ia1) + Nwann
iend1= orbital_start(ia1+1)- 1 + Nwann
istart2= orbital_start(ia2) + Nwann
iend2= orbital_start(ia2+1)- 1 + Nwann
H01 ( istart1:iend1, istart2:iend2) &
= H01 ( istart1:iend1, istart2:iend2) &
+ HmnR( istart1:iend1, istart2:iend2, iR)*ratio/ndegen(iR)* fac
endif ! soc
enddo ! ia2
enddo ! ia1
enddo ! iR
do coo2=1,Num_wann!*(i1-1),Num_wann*i1
do coo1=1,Num_wann!*(i2-1),Num_wann*i2
tmp=H01(coo1,coo2)
if(abs(tmp)/eV2Hartree >eps6) then
!> only store the upper triangle matrix for mkl_dss purpose
j_t=coo1+(i1-1)*Num_wann
i_t=coo2+(i2-1)*Num_wann
!if (i_t>=j_t) then
ncoo=ncoo+1
jcoo(ncoo)=j_t
icoo(ncoo)=i_t
acoo(ncoo)=tmp*exp(zi*k(3))
!endif
!> add the transpose-conjugate part
i_t=coo1+(i1-1)*Num_wann
j_t=coo2+(i2-1)*Num_wann
!if (i_t>=j_t) then
ncoo=ncoo+1
jcoo(ncoo)=j_t
icoo(ncoo)=i_t
acoo(ncoo)=conjg(tmp)*exp(-zi*k(3)) !0.01
!endif
end if
end do
end do
enddo ! i2
enddo ! i1
if (cpuid.eq.0) write(stdout,*) 'nnz, nnzmax after H01:', ncoo, nnz
if (nnz< ncoo) then
write(*, *)'Please increase nnzmax in the sparse.f90'
write(*,*) 'nnz, nnz after H01 : ', ncoo, nnz
stop
endif
nnz= ncoo
return
end subroutine ham_3Dlandau_sparse1
subroutine ham_3Dlandau_sparseHR(nnz, Ndimq, Nq, k, acoo,jcoo,icoo)
!> We generate the TB hamiltonian under magnetic field with sparse formated HmnR.
!> At the begining, we generate coo-formated hamlitonian. However, in the end,
!> the coo-formated hamlitonian will be transform into csr format
use para
implicit none
!> input: nnz is the maximum number of non-zeros entries
!> output: nnz is the number of non-zeros entries of acoo
integer, intent(inout) :: nnz
integer, intent(in) :: Ndimq
integer, intent(in) :: Nq
real(dp), intent(in) :: k(3)
!> output hamiltonian stored as COO sparse matrix format
complex(dp), intent(inout) :: acoo(nnz)
integer, intent(inout) :: jcoo(nnz)
integer, intent(inout) :: icoo(nnz)
integer,allocatable :: rxyz(:,:),nonzero_counter(:,:,:)
! loop index
integer :: i1, i2,i, i_t, j_t,i3
integer :: coo1,coo2,ncoo
! loop index
integer :: iR,ims
! index used to sign irvec
real(dp) :: ia,ib,ic
integer :: ia1, ia2
integer :: inew_ic,inew_ia,inew_ib
!> new index used to sign irvec
real(dp) :: new_ia,new_ib,new_ic
!> wave vector k times lattice vector R
real(Dp) :: kdotr
real(dp) :: phase
complex(dp) :: ratio
complex(dp) :: fac,tmp
real(dp) :: Rp1(3), Rp2(3), R1(3), R2(3)
real(dp) :: Ri(3), Rj(3), tau1(3), tau2(3)
!> calculate phase induced by magnetic field
real(dp),external :: phase1,phase2
!> a work array for coocsr
integer(li), allocatable, target :: iwk (:)
if(export_maghr) then
allocate(rxyz(nnz,3))
allocate(nonzero_counter(-ijmax:ijmax,-ijmax:ijmax,-ijmax:ijmax))
nonzero_counter=0
endif
!> calculate intra-hopping
acoo=zzero
ncoo=0
tmp=0d0
! i1 column index, sweep magnetic supercells
do i1=1, Nq
! i2 row index, sweep magnetic supercells
do i2=1, Nq
if (abs(i2-i1)> ijmax) cycle
!> sum over R points to get H(k1, k2)
do ims=1,splen
ia=dble(hirv(1,ims))
ib=dble(hirv(2,ims))
ic=dble(hirv(3,ims))
!> new lattice
call latticetransform(ia, ib, ic, new_ia, new_ib, new_ic)
inew_ic= int(new_ic)
inew_ib=int(new_ib);inew_ia=int(new_ia)
if (inew_ic /= (i2-i1)) cycle
!> exp(i k.R)
kdotr= k(1)*new_ia+ k(2)*new_ib
ratio= cos(2d0*pi*kdotr)+ zi*sin(2d0*pi*kdotr)
Rp1= (i1-1)*Ruc_new
Rp2= new_ia*Rua_new+ new_ib*Rub_new+ (i2-1)*Ruc_new
ia1=hicoo(ims)
ia2=hjcoo(ims)
R1= Origin_cell%wannier_centers_cart(:, ia1)
R2= Origin_cell%wannier_centers_cart(:, ia2)
call rotate(R1, tau1)
call rotate(R2, tau2)
Ri= Rp1+ tau1
Rj= Rp2+ tau2
phase=phase2(Ri,Rj)
fac= cos(phase)+ zi*sin(phase)
tmp=hacoo(ims)*ratio/ndegen(iR)* fac
if(abs(tmp)/eV2Hartree > 1e-6) then
ncoo=ncoo+1
icoo(ncoo)=hicoo(ims)+(i1-1)*Num_wann
jcoo(ncoo)=hjcoo(ims)+(i2-1)*Num_wann
acoo(ncoo)=acoo(ncoo)+tmp
!if (export_maghr) then
! rxyz(ncoo,:)=[inew_ia,inew_ib,inew_ic]
! nonzero_counter(inew_ia,inew_ib,inew_ic)=1
!endif
endif
enddo ! iR
enddo ! i2
enddo ! i1
!> if we want to calculate the Chern number, we use open boundary such
!> that we can get the Chern number from the edge states.
if(landau_chern_calc) goto 239
do i1=1, Nq
! i2 row index, sweep magnetic supercells
do i2=1, Nq
if (abs(i2-i1+Nq)> ijmax) cycle
!> sum over R points to get H(k1, k2)
do ims=1,splen
ia=dble(hirv(1,ims))
ib=dble(hirv(2,ims))
ic=dble(hirv(3,ims))
!> new lattice
call latticetransform(ia, ib, ic, new_ia, new_ib, new_ic)
inew_ic= int(new_ic)
inew_ib=int(new_ib);inew_ia=int(new_ia)
if (inew_ic /= (i2-i1+Nq)) cycle
!> exp(i k.R)
kdotr= k(1)*new_ia+ k(2)*new_ib+ k(3)
ratio= cos(2d0*pi*kdotr)+ zi*sin(2d0*pi*kdotr)
Rp1= (i1-1)*Ruc_new
Rp2= new_ia*Rua_new+ new_ib*Rub_new+ (i2-1+Nq)*Ruc_new
ia1=hicoo(ims)
ia2=hjcoo(ims)
R1= Origin_cell%wannier_centers_cart(:, ia1)
R2= Origin_cell%wannier_centers_cart(:, ia2)
call rotate(R1, tau1)
call rotate(R2, tau2)
Ri= Rp1+ tau1
Rj= Rp2+ tau2
phase=phase2(Ri,Rj)
fac= cos(phase)+ zi*sin(phase)
tmp=hacoo(ims)*ratio/ndegen(iR)* fac
if(abs(tmp)/eV2Hartree > 1e-6) then
ncoo=ncoo+1
icoo(ncoo)=hicoo(ims)+(i1-1)*Num_wann
jcoo(ncoo)=hjcoo(ims)+(i2-1)*Num_wann
acoo(ncoo)=acoo(ncoo)+tmp
if(export_maghr) then
rxyz(ncoo,:)=[inew_ia,inew_ib,inew_ic]
nonzero_counter(inew_ia,inew_ib,inew_ic)=1
endif
endif
tmp=conjg(tmp)
if(abs(tmp)/eV2Hartree > 1e-6) then
ncoo=ncoo+1
icoo(ncoo)=hjcoo(ims)+(i2-1)*Num_wann
jcoo(ncoo)=hicoo(ims)+(i1-1)*Num_wann
acoo(ncoo)=acoo(ncoo)+tmp
!if (export_maghr) then
! rxyz(ncoo,:)=-[inew_ia,inew_ib,inew_ic]
! nonzero_counter(-inew_ia,-inew_ib,-inew_ic)=1
!endif
endif
enddo ! iR
enddo ! i2
enddo ! i1
239 continue
if (cpuid.eq.0) write(stdout,*) 'nnz, nnzmax after H01:', ncoo, nnz
if (nnz< ncoo) then
write(*, *)'>>> Error : Please increase nnz in the sparse.f90'
write(*, *) 'nnz, nnzmax after H01 : ', ncoo, nnz
stop
endif
nnz= ncoo
!> usually, we don't export magnetic hr file
if(export_maghr) then
outfileindex= outfileindex+ 1
open(unit=outfileindex,file='hr_mag.dat')
write(outfileindex,*) '!magnetic supercell hr automated'
write(outfileindex,*) ncoo
write(outfileindex,*) ndimq
write(outfileindex,*) sum(nonzero_counter)
write(outfileindex,*) ('1 ',i=1,sum(nonzero_counter))
do i1=-ijmax,ijmax
do i2=-ijmax,ijmax
do i3=-ijmax,ijmax
if(nonzero_counter(i1,i2,i3)==0) cycle
do i=1,ncoo
if(sum(abs(rxyz(i,:)-[i1,i2,i3]))/=0) cycle
write(outfileindex,*) rxyz(i,:),icoo(i),jcoo(i),real(acoo(i)),imag(acoo(i))
end do
end do
end do
end do
endif
if(export_maghr) deallocate(rxyz)
return
end subroutine ham_3Dlandau_sparseHR
!use this to make eigs from coo
!> use solver from ARPACK
subroutine sparse_landau_level_B
use para
use sparse
implicit none
!> magnetic supercell size, perpendicular to the magnetic field
!> Ndimq= Nq* Num_wann
integer :: Nq, Ndimq, Nmag, Nmag1
integer :: ia1, ia2, ib, i, j, kk1, kk2, ie, iq, ierr, ig
real(dp) :: B0, B0Tesla, B0Tesla_quantumflux_magsupcell, theta, dis, dis1
! wave vector
real(dp) :: k3(3)
!> dim= Ndimq, knv3
real(dp), allocatable :: mag(:), mag_Tesla(:), W(:), eigv(:, :), eigv_mpi(:, :)
real(dp) :: time_start, time_end, time_start0
!> dim= Ndimq*Ndimq
!> maximum number of non-zeros matrix elements
integer :: nnzmax, nnz
complex(dp), allocatable :: acoo(:)
integer, allocatable :: icoo(:), jcoo(:)
!number of ARPACK eigenvalues to be obtained
integer :: neval
! number of Arnoldi vectors
integer :: nvecs
!shift-invert sigma
complex(dp) :: sigma
!> time measurement
real(dp) :: time1, time2, time3, times, timee
!> eigenvector of the sparse matrix acoo. Dim=(Ndimq, neval)
complex(dp), allocatable :: psi(:)
!> the dimension of zeigv is (Ndimq, nvecs)
!> however, only (1:Ndimq, 1:neval) are useful
complex(dp), allocatable :: zeigv(:, :)
!> print the weight for the Selected_WannierOrbitals
real(dp), allocatable :: dos_selected(:, :, :), dos_selected_mpi(:, :, :)
Nq= Magq
Ndimq= Num_wann* Nq
if (NumSelectedEigenVals==0) NumSelectedEigenVals=OmegaNum
neval=NumSelectedEigenVals
if (neval>=Ndimq) neval= Ndimq- 2
!> ncv
nvecs=int(2*neval)
if (nvecs<50) nvecs= 50
if (nvecs>Ndimq) nvecs= Ndimq
sigma=(1d0,0d0)*E_arc
Nmag= Magp-1
Nmag1=Nmag/1
nnzmax= Num_wann*(2*ijmax+1)*Ndimq+Ndimq
if(Is_Sparse_Hr) nnzmax=splen*nq+Ndimq
allocate( acoo(nnzmax), stat=ierr)
call printallocationinfo('acoo', ierr)
allocate( jcoo(nnzmax), stat=ierr)
call printallocationinfo('jcoo', ierr)
allocate( icoo(nnzmax), stat=ierr)
call printallocationinfo('icoo', ierr)
allocate( W( neval))
allocate( eigv( neval, Nmag1+1), stat=ierr)
call printallocationinfo('eigv', ierr)
allocate( eigv_mpi( neval, Nmag1+1), stat=ierr)
call printallocationinfo('eigv_mpi', ierr)
allocate( mag(Nmag+1), mag_Tesla(Nmag+ 1))
allocate( psi(ndimq))
allocate( zeigv(ndimq,nvecs), stat=ierr)
call printallocationinfo('zeigv', ierr)
allocate( dos_selected (neval, Nmag+ 1, NumberofSelectedOrbitals_groups), stat=ierr)
call printallocationinfo('dos_selected', ierr)
allocate( dos_selected_mpi (neval, Nmag+ 1, NumberofSelectedOrbitals_groups), stat=ierr)
call printallocationinfo('dos_selected_mpi', ierr)
zeigv= 0d0
dos_selected= 0d0
dos_selected_mpi= 0d0
mag= 0d0; mag_Tesla= 0d0
eigv_mpi= 0d0
eigv = 0d0
acoo=0d0
jcoo=0
icoo=0
if (cpuid.eq.0) write(stdout,*) 'sigma=',sigma,Ndimq,NumSelectedEigenVals
!> deal with the magnetic field
!> first transform the Bx By into B*Cos\theta, B*Sin\theta
if (abs(By)<1e-8) then
if (Bx<0) then
theta= pi
else
theta= 0d0
endif
elseif (By>0) then
theta = atan(Bx/By)
else
theta = atan(Bx/By)+ pi
endif
!> the magnetic flux in the magnetic supercell is 2*pi*Magp
B0= 2d0*pi/dble(Nq)
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/1.6021766208*1E19/MagneticSuperProjectedArea*1E20
B0Tesla= 6.62607004*1E-34/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
k3=Single_KPOINT_3D_DIRECT
!> calculate the landau levels along special k line
!> The flux in the unit cell changes from 0 to 2*pi
do ib=1, Nmag+ 1
mag(ib)= B0* (ib-1)
mag_Tesla(ib)= B0Tesla_quantumflux_magsupcell* (ib-1)
enddo
time_start= 0d0
time_start0= 0d0
call now(time_start0)
time_start= time_start0
time_end = time_start0
do ib=1+ cpuid, Nmag1+1, num_cpu
call now(times)
!Bx= mag(ib)* Cos(theta)
!By= mag(ib)* Sin(theta)
Bx= -mag(ib)
By= 0d0
nnz= nnzmax
if (cpuid.eq.0) &
write(stdout, '(2a, i9, " /", i10, a, f10.1, "s", a, f10.1, "s")') &
'In sparse_landau_level_B ', ' ib/Nmag ', ib, Nmag+1, &
' time elapsed: ', time_end-time_start0, &
' time left: ', ((Nmag-ib)/num_cpu)*(time_end-time_start)
call now(time_start)
call now(time1)
if(Is_Sparse_Hr) then
call ham_3Dlandau_sparseHR(nnz, Ndimq, Nq, k3, acoo,jcoo,icoo)
else
call ham_3Dlandau_sparse1(nnz, Ndimq, Nq, k3, acoo,jcoo,icoo)
endif
acoo=acoo/eV2Hartree
call now(time2)
!> diagonalization by call zheev in lapack
W= 0d0
#if defined (INTELMKL)
call arpack_sparse_coo_eigs(Ndimq,nnzmax,nnz,acoo,jcoo,icoo,neval,nvecs,W,sigma, zeigv, LandauLevel_wavefunction_calc)
#endif
call now(time3)
eigv(:, ib)= W
if (cpuid==0)write(stdout, '(a, f20.2, a)')' >> Time cost for constructing H: ', time2-time1, ' s'
if (cpuid==0)write(stdout, '(a, f20.2, a)')' >> Time cost for diagonalize H: ', time3-time2, ' s'
if (LandauLevel_wavefunction_calc) then
do ie= 1, neval
psi(:)= zeigv(:, ie)
do ig=1, NumberofSelectedOrbitals_groups
do iq=1, Nq
do i= 1, NumberofSelectedOrbitals(ig)
j= Num_wann*(iq-1)+ Selected_WannierOrbitals(ig)%iarray(i)
dos_selected(ie, ib, ig)= dos_selected(ie, ib, ig)+ abs(psi(j))**2
enddo ! sweep the selected orbitals
enddo ! iq sweep the magnetic supercell
enddo ! ig group of selected orbitals
enddo ! ib sweep the eigenvalue
endif
call now(time_end)
enddo !ib
#if defined (MPI)
call mpi_allreduce(eigv,eigv_mpi,size(eigv),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
call mpi_allreduce(dos_selected, dos_selected_mpi,size(dos_selected),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
#else
eigv_mpi= eigv
dos_selected_mpi= dos_selected
#endif
outfileindex= outfileindex+ 1
if (cpuid.eq.0) then
open(unit=outfileindex, file='landaulevel_B.dat')
if (LandauLevel_wavefunction_calc) then
write(outfileindex, '("#", a14, 2a15, a)')'Phi per cell', 'B (Tesla)', ' Eig of LL', ' Weight on the selected orbitals'
do j=1, neval
do i=1, Nmag1+1
write(outfileindex,'(4f19.7)')mag(i)/2d0/pi, mag_Tesla(i), eigv_mpi(j, i), &
(dos_selected_mpi(j, i, ig), ig=1, NumberofSelectedOrbitals_groups)
enddo
write(outfileindex , *)''
enddo
else
write(outfileindex, '("#", a14, 2a15)')'Phi per cell', 'B (Tesla)', ' Eig of LL'
do j=1, neval
do i=1, Nmag1+1
write(outfileindex,'(3f19.7)')mag(i)/2d0/pi, mag_Tesla(i), eigv_mpi(j, i)
enddo
write(outfileindex , *)''
enddo
endif
close(outfileindex)
write(stdout,*) 'calculate landau level done'
endif
outfileindex= outfileindex+ 1
if (cpuid==0) then
open(unit=outfileindex, file='landaulevel_B.gnu')
write(outfileindex, '(a)')"set encoding iso_8859_1"
write(outfileindex, '(a)')'#set terminal postscript enhanced color'
write(outfileindex, '(a)')"#set output 'landaulevel_B.eps'"
write(outfileindex, '(3a)')'#set terminal pngcairo truecolor enhanced', &
' font ",60" size 1920, 1680'
write(outfileindex, '(3a)')'set terminal png truecolor enhanced', &
' font ",60" size 1920, 1680'
write(outfileindex, '(a)')"set output 'landaulevel_B.png'"
write(outfileindex, '(a)')'unset key'
write(outfileindex, '(a)')'set pointsize 0.8'
write(outfileindex, '(a)')'set palette rgb 21,22,23'
write(outfileindex, '(a)')'set border lw 3 '
write(outfileindex, '(a)')'#set xtics font ",36"'
write(outfileindex, '(a)')'#set ytics font ",36"'
write(outfileindex, '(a)')'#set xlabel font ",36"'
write(outfileindex, '(a)')'set xlabel "Phi per cell"'
write(outfileindex,*) '#set xlabel "{/Symbol F}/{/Symbol F}_0"'
write(outfileindex, '(a)')'set ylabel "Energy (eV)"'
write(outfileindex, '(a, i6,a)')'set title "Landau level with Nq=', Nq, '"'
write(outfileindex, '(a)')'#set xtics offset 0, -1'
write(outfileindex, '(a)')'#set ylabel offset -1, 0 '
write(outfileindex, '(a)')'rgb(r,g,b) = int(r)*65536 + int(g)*256 + int(b)'
if (LandauLevel_wavefunction_calc) then
write(outfileindex, '(2a)')"plot 'landaulevel_B.dat' u 1:3:(rgb(255,255-255*$3, 3)) ", &
"w p pt 7 ps 1 lc rgb variable"
else
write(outfileindex, '(2a)')"plot 'landaulevel_B.dat' u 1:3", &
" w p pt 7 ps 1"
endif
close(outfileindex)
endif
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
deallocate( acoo)
deallocate( jcoo)
deallocate( icoo)
deallocate( W)
deallocate( eigv)
deallocate( eigv_mpi)
deallocate( mag)
return
end subroutine sparse_landau_level_B
subroutine sparse_landau_level_k
use para
use sparse
implicit none
!> magnetic supercell size, perpendicular to the magnetic field
!> Ndimq= Nq* Num_wann
integer :: Nq, Ndimq
integer :: Nmag,Nmag1
!> some temporary integers
integer :: ik, ia1, ia2, i, j, ierr, ib, iq, ig
real(dp) :: B0, B0Tesla, B0Tesla_quantumflux_magsupcell
real(dp) :: theta
real(dp) :: t1, temp, dis, dis1
! wave vector
real(dp) :: k3(3)
!> dim= Ndimq, knv3
real(dp), allocatable :: W(:)
real(dp), allocatable :: eigv(:, :)
real(dp), allocatable :: eigv_mpi(:, :)
real(dp) :: emin, emax
real(dp) :: time_start, time_end, time_start0
!> dim= Ndimq*Ndimq
integer :: nnzmax, nnz
complex(dp), allocatable :: acoo(:)
integer, allocatable :: jcoo(:)
integer, allocatable :: icoo(:)
!> eigenvector of the sparse matrix acoo. Dim=(Ndimq, neval)
complex(dp), allocatable :: psi(:)
complex(dp), allocatable :: zeigv(:, :)
!> print the weight for the Selected_WannierOrbitals
real(dp), allocatable :: dos_selected(:, :, :)
real(dp), allocatable :: dos_selected_mpi(:, :, :)
real(dp), allocatable :: dos_l_selected(:, :, :)
real(dp), allocatable :: dos_l_selected_mpi(:, :, :)
real(dp), allocatable :: dos_r_selected(:, :, :)
real(dp), allocatable :: dos_r_selected_mpi(:, :, :)
!number of ARPACK eigenvalues
integer :: neval
! number of Arnoldi vectors
integer :: nvecs
!shift-invert sigma
complex(dp) :: sigma
!> time measurement
real(dp) :: time1, time2, time3
Nq= Magq
Nmag= Nq
Ndimq= Num_wann* Nq
nnzmax= Num_wann*(2*ijmax+1)*Ndimq+Ndimq
if(Is_Sparse_Hr) nnzmax=splen*nq+Ndimq
if (NumSelectedEigenVals==0) NumSelectedEigenVals=OmegaNum
neval=NumSelectedEigenVals
if (neval>=Ndimq) neval= Ndimq- 2
!> ncv
nvecs=int(2*neval)
if (nvecs<50) nvecs= 50
if (nvecs>Ndimq) nvecs= Ndimq
sigma=(1d0,0d0)*E_arc
allocate( acoo(nnzmax), stat= ierr)
call printallocationinfo('acoo', ierr)
allocate( jcoo(nnzmax), stat= ierr)
call printallocationinfo('jcoo', ierr)
allocate( icoo(nnzmax), stat= ierr)
call printallocationinfo('icoo', ierr)
allocate( W( neval), stat= ierr)
allocate( eigv( neval, nk3_band))
call printallocationinfo('eigv', ierr)
allocate( eigv_mpi( neval, nk3_band), stat= ierr)
call printallocationinfo('eigv_mpi', ierr)
allocate( psi(ndimq))
allocate( zeigv(ndimq,nvecs), stat= ierr)
call printallocationinfo('zeigv', ierr)
allocate( dos_selected (neval, nk3_band, NumberofSelectedOrbitals_groups), stat= ierr)
call printallocationinfo('dos_selected', ierr)
allocate( dos_selected_mpi (neval, nk3_band, NumberofSelectedOrbitals_groups), stat= ierr)
call printallocationinfo('dos_selected_mpi', ierr)
if (landau_chern_calc) then
allocate( dos_l_selected (neval, nk3_band, NumberofSelectedOrbitals_groups), stat= ierr)
call printallocationinfo('dos_l_selected', ierr)
allocate( dos_l_selected_mpi (neval, nk3_band, NumberofSelectedOrbitals_groups), stat= ierr)
call printallocationinfo('dos_l_selected_mpi', ierr)
allocate( dos_r_selected (neval, nk3_band, NumberofSelectedOrbitals_groups), stat= ierr)
call printallocationinfo('dos_r_selected', ierr)
allocate( dos_r_selected_mpi (neval, nk3_band, NumberofSelectedOrbitals_groups), stat= ierr)
call printallocationinfo('dos_r_selected_mpi', ierr)
dos_l_selected= 0d0
dos_l_selected_mpi= 0d0
dos_r_selected= 0d0
dos_r_selected_mpi= 0d0
endif
dos_selected= 0d0
dos_selected_mpi= 0d0
if (cpuid==0) then
write(stdout, '(a, I16)')' nnzmax is ', &
nnzmax
write(stdout, '(a, I16)')' Dimension of magnetic supercell Ndimq is ', &
ndimq
write(stdout, '(a, I16)')' Number of ritz vectors nvecs is ', &
nvecs
write(stdout, '(a, f20.1, a)')' Memory for acoo, icoo, jcoo is ', &
nnzmax*(16d0+4d0+4d0)/1024d0/1024d0, ' MB'
write(stdout, '(a, f20.1, a)')' Memory for zeigv matrix is ', &
ndimq/1024d0/1024d0*nvecs*16d0, ' MB'
write(stdout, '(a, f20.1, a)')' Memory for eigv matrix is ', &
neval*nk3_band*16d0/1024d0/1024d0, ' MB'
if (landau_chern_calc) then
write(stdout, '(a, f20.1, a)')' Memory for dos_selected is ', &
neval*nk3_band*6*NumberofSelectedOrbitals_groups*8d0/1024d0/1024d0, ' MB'
else
write(stdout, '(a, f20.1, a)')' Memory for dos_selected is ', &
neval*nk3_band*2*NumberofSelectedOrbitals_groups*8d0/1024d0/1024d0, ' MB'
endif
write(stdout, *)' '