-
Notifications
You must be signed in to change notification settings - Fork 0
/
berrycurvature.f90
executable file
·2497 lines (2076 loc) · 91.5 KB
/
berrycurvature.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
subroutine berry_curvarture_singlek_EF(k, mu, Omega_x, Omega_y, Omega_z)
!> Calculate Berry curvature for a sigle k point
!> The Fermi distribution is determined by the Fermi level E_arc
!> ref : Physical Review B 74, 195118(2006)
!
!> eqn (34)
!
!> Jun. 11 2018 by Quansheng Wu @ ETHZ
!> Try to write some simple subroutines.
!> on the output, only the real part is meaningfull.
! Copyright (c) 2018 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
!> input parameter, k is in unit of reciprocal lattice vectors
!> mu: chemical potential which can be tuned by the back gate
real(dp), intent(in) :: k(3), mu
complex(dp), intent(out) :: Omega_x(Num_wann)
complex(dp), intent(out) :: Omega_y(Num_wann)
complex(dp), intent(out) :: Omega_z(Num_wann)
integer :: m, n, i
real(dp), allocatable :: W(:)
real(dp) :: beta_fake
complex(dp), allocatable :: Amat(:, :), DHDk(:, :, :), DHDkdag(:, :, :)
complex(dp), allocatable :: UU(:, :), UU_dag(:, :), Hamk_bulk(:, :)
complex(dp), allocatable :: velocity_wann(:, :, :), velocity_Ham(:, :, :)
!> Fermi-Dirac distribution
real(dp), external :: fermi
allocate(W(Num_wann))
allocate(velocity_wann(Num_wann, Num_wann, 3), velocity_Ham(Num_wann, Num_wann, 3))
allocate(UU(Num_wann, Num_wann), UU_dag(Num_wann, Num_wann), Hamk_bulk(Num_wann, Num_wann))
allocate(Amat(Num_wann, Num_wann), DHDk(Num_wann, Num_wann, 3), DHDkdag(Num_wann, Num_wann, 3))
W=0d0; velocity_wann= 0d0; UU= 0d0; UU_dag= 0d0; velocity_Ham= 0d0
! calculation bulk hamiltonian by a direct Fourier transformation of HmnR
call ham_bulk_atomicgauge(k, Hamk_bulk)
!> diagonalization by call zheev in lapack
UU=Hamk_bulk
call eigensystem_c( 'V', 'U', Num_wann, UU, W)
!call zhpevx_pack(hamk_bulk,Num_wann, W, UU)
UU_dag= conjg(transpose(UU))
call dHdk_atomicgauge(k, velocity_wann)
!> unitility rotate velocity
do i=1, 3
call mat_mul(Num_wann, velocity_wann(:, :, i), UU, Amat)
call mat_mul(Num_wann, UU_dag, Amat, velocity_Ham(:, :, i))
enddo
Omega_x=0d0;Omega_y=0d0; Omega_z=0d0
do m= 1, Num_wann
do n= 1, Num_wann
if (abs(W(m)-W(n))<eps6) cycle
Omega_x(m)= Omega_x(m)+ velocity_Ham(n, m, 2)*velocity_Ham(m, n, 3)/((W(m)-W(n))**2)
Omega_y(m)= Omega_y(m)+ velocity_Ham(n, m, 3)*velocity_Ham(m, n, 1)/((W(m)-W(n))**2)
Omega_z(m)= Omega_z(m)+ velocity_Ham(n, m, 1)*velocity_Ham(m, n, 2)/((W(m)-W(n))**2)
enddo ! m
enddo ! n
Omega_x= -Omega_x*2d0*zi
Omega_y= -Omega_y*2d0*zi
Omega_z= -Omega_z*2d0*zi
!> consider the Fermi-distribution according to the broadening Earc_eta
if (Eta_Arc<eps6) Eta_Arc= eps6
beta_fake= 1d0/Eta_Arc
do m= 1, Num_wann
Omega_x(m)= Omega_x(m)*fermi(W(m)-mu, beta_fake)
Omega_y(m)= Omega_y(m)*fermi(W(m)-mu, beta_fake)
Omega_z(m)= Omega_z(m)*fermi(W(m)-mu, Beta_fake)
enddo
return
end subroutine berry_curvarture_singlek_EF
subroutine berry_curvarture_singlek_numoccupied_slab_total(k, Omega_z)
!> Calculate Berry curvature for a sigle k point
!> The Fermi distribution is determined by the NumOccupied bands, not the Fermi level
!> ref : Physical Review B 74, 195118(2006)
!
!> eqn (11), we only calculate xy
!
!> Aug. 06 2018 by Quansheng Wu
!> Try to write some simple subroutines.
! Copyright (c) 2018 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
!> input parameter, k is in unit of reciprocal lattice vectors
real(dp), intent(in) :: k(2)
complex(dp), intent(out) :: Omega_z(1)
integer :: m, n, Mdim, i1, i2
real(dp), allocatable :: W(:)
complex(dp), allocatable :: Amat(:, :), DHDk(:, :, :), DHDkdag(:, :, :)
complex(dp), allocatable :: UU(:, :), UU_dag(:, :), Hamk_slab(:, :)
complex(dp), allocatable :: vx(:, :), vy(:, :)
complex(dp), allocatable :: Vij_x(:, :, :), Vij_y(:, :, :)
!> leading order of the hamiltonian matrix for slab system specified with Nslab
Mdim = Num_wann*Nslab
allocate(W(Mdim))
allocate(vx(Mdim, Mdim), vy(Mdim, Mdim))
allocate(UU(Mdim, Mdim), UU_dag(Mdim, Mdim), Hamk_slab(Mdim, Mdim))
allocate(Amat(Mdim, Mdim), DHDk(Mdim, Mdim, 3), DHDkdag(Mdim, Mdim, 3))
allocate(Vij_x(-ijmax:ijmax, Num_wann, Num_wann))
allocate(Vij_y(-ijmax:ijmax, Num_wann, Num_wann))
W=0d0; vx= 0d0; vy= 0d0; UU= 0d0; UU_dag= 0d0
call ham_qlayer2qlayer_velocity(k, Vij_x, Vij_y)
do i1=1, nslab
! i2 row index
do i2=1, nslab
if (abs(i2-i1).le.ijmax)then
vx((i2-1)*Num_wann+1:(i2-1)*Num_wann+Num_wann,&
(i1-1)*Num_wann+1:(i1-1)*Num_wann+Num_wann )&
= Vij_x(i2-i1,1:Num_wann,1:Num_wann)
vy((i2-1)*Num_wann+1:(i2-1)*Num_wann+Num_wann,&
(i1-1)*Num_wann+1:(i1-1)*Num_wann+Num_wann )&
= Vij_y(i2-i1,1:Num_wann,1:Num_wann)
endif
enddo ! i2
enddo ! i1
! calculation slab hamiltonian by a direct Fourier transformation of HmnR
call ham_slab(k, Hamk_slab)
!> diagonalization by call zheev in lapack
UU=Hamk_slab
call eigensystem_c( 'V', 'U', Mdim, UU, W)
!call zhpevx_pack(hamk_slab,Mdim, W, UU)
UU_dag= conjg(transpose(UU))
!> unitility rotate velocity
call mat_mul(Mdim, vx, UU, Amat)
call mat_mul(Mdim, UU_dag, Amat, vx)
call mat_mul(Mdim, vy, UU, Amat)
call mat_mul(Mdim, UU_dag, Amat, vy)
Omega_z=0d0
do m= 1, NumOccupied*Nslab
do n= NumOccupied*Nslab+1, Mdim
Omega_z(1)= Omega_z(1)+ vx(n, m)*vy(m, n)/((W(m)-W(n))**2)
enddo ! m
enddo ! n
Omega_z= -aimag(Omega_z*2d0)/Angstrom2atomic**2
return
end subroutine berry_curvarture_singlek_numoccupied_slab_total
subroutine berry_curvarture_singlek_numoccupied_total(k, Omega_x, Omega_y, Omega_z)
!> Calculate Berry curvature for a sigle k point using Kubo-formula
!> The Fermi distribution is determined by the NumOccupied bands, not the Fermi level
!> ref : Physical Review B 74, 195118(2006)
!
!> eqn (11), we only calculate yz, zx, xy
!
!> Jun. 25 2018 by Quansheng Wu @ airplane CA781 from Beijing to Zurich
!> Try to write some simple subroutines.
! Copyright (c) 2018 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
!> input parameter, k is in unit of reciprocal lattice vectors
real(dp), intent(in) :: k(3)
complex(dp), intent(out) :: Omega_x
complex(dp), intent(out) :: Omega_y
complex(dp), intent(out) :: Omega_z
integer :: m, n, i
real(dp), allocatable :: W(:)
complex(dp), allocatable :: Amat(:, :)
complex(dp), allocatable :: UU(:, :), UU_dag(:, :), Hamk_bulk(:, :)
complex(dp), allocatable :: velocity_wann(:, :, :), velocity_Ham(:, :, :)
allocate(W(Num_wann))
allocate(velocity_wann(Num_wann, Num_wann, 3), velocity_Ham(Num_wann, Num_wann, 3))
allocate(UU(Num_wann, Num_wann), UU_dag(Num_wann, Num_wann), Hamk_bulk(Num_wann, Num_wann))
allocate(Amat(Num_wann, Num_wann))
W=0d0; velocity_wann= 0d0; UU= 0d0; UU_dag= 0d0; velocity_Ham= 0d0
! calculation bulk hamiltonian by a direct Fourier transformation of HmnR
!call ham_bulk_latticegauge(k, Hamk_bulk)
call ham_bulk_atomicgauge (k, Hamk_bulk)
!> diagonalization by call zheev in lapack
UU=Hamk_bulk
call eigensystem_c( 'V', 'U', Num_wann, UU, W)
!call zhpevx_pack(hamk_bulk,Num_wann, W, UU)
UU_dag= conjg(transpose(UU))
call dHdk_atomicgauge(k, velocity_wann)
!> unitility rotate velocity
do i=1, 3
call mat_mul(Num_wann, velocity_wann(:, :, i), UU, Amat)
call mat_mul(Num_wann, UU_dag, Amat, velocity_Ham(:, :, i))
enddo
Omega_x=0d0; Omega_y=0d0; Omega_z=0d0
do m= 1, NumOccupied !> sum over valence band
do n= NumOccupied+1, Num_wann !> sum over conduction band
Omega_x= Omega_x+ velocity_Ham(n, m, 2)*velocity_Ham(m, n, 3)/((W(m)-W(n))**2)
Omega_y= Omega_y+ velocity_Ham(n, m, 3)*velocity_Ham(m, n, 1)/((W(m)-W(n))**2)
Omega_z= Omega_z+ velocity_Ham(n, m, 1)*velocity_Ham(m, n, 2)/((W(m)-W(n))**2)
enddo ! m
enddo ! n
Omega_x= aimag(Omega_x*2d0)
Omega_y= aimag(Omega_y*2d0)
Omega_z= aimag(Omega_z*2d0)
deallocate(W, UU, UU_dag, hamk_bulk)
deallocate(Amat, velocity_wann, velocity_Ham)
return
end subroutine berry_curvarture_singlek_numoccupied_total
subroutine berry_curvarture_singlek_numoccupied(k, Omega_x, Omega_y, Omega_z)
!> Calculate Berry curvature for a sigle k point
!> The Fermi distribution is determined by the NumOccupied bands, not the Fermi level
!> ref : Physical Review B 74, 195118(2006)
!
!> eqn (30), we only calculate yz, zx, xy
!
!> Jun. 25 2018 by Quansheng Wu @ airplane CA781 from Beijing to Zurich
!> Try to write some simple subroutines.
! Copyright (c) 2018 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
!> input parameter, k is in unit of reciprocal lattice vectors
real(dp), intent(in) :: k(3)
complex(dp), intent(out) :: Omega_x(Num_wann)
complex(dp), intent(out) :: Omega_y(Num_wann)
complex(dp), intent(out) :: Omega_z(Num_wann)
integer :: m, n, i
real(dp), allocatable :: W(:)
complex(dp), allocatable :: Amat(:, :), DHDk(:, :, :), DHDkdag(:, :, :)
complex(dp), allocatable :: UU(:, :), UU_dag(:, :), Hamk_bulk(:, :)
complex(dp), allocatable :: velocity_wann(:, :, :), velocity_Ham(:, :, :)
allocate(W(Num_wann))
allocate(velocity_wann(Num_wann, Num_wann, 3), velocity_Ham(Num_wann, Num_wann, 3))
allocate(UU(Num_wann, Num_wann), UU_dag(Num_wann, Num_wann), Hamk_bulk(Num_wann, Num_wann))
allocate(Amat(Num_wann, Num_wann), DHDk(Num_wann, Num_wann, 3), DHDkdag(Num_wann, Num_wann, 3))
W=0d0; velocity_wann= 0d0; UU= 0d0; UU_dag= 0d0; velocity_Ham= 0d0
! calculation bulk hamiltonian by a direct Fourier transformation of HmnR
call ham_bulk_atomicgauge (k, Hamk_bulk)
!> diagonalization by call zheev in lapack
UU=Hamk_bulk
call eigensystem_c( 'V', 'U', Num_wann, UU, W)
!call zhpevx_pack(hamk_bulk,Num_wann, W, UU)
UU_dag= conjg(transpose(UU))
call dHdk_atomicgauge(k, velocity_wann)
!> unitility rotate velocity
do i=1, 3
call mat_mul(Num_wann, velocity_wann(:, :, i), UU, Amat)
call mat_mul(Num_wann, UU_dag, Amat, velocity_Ham(:, :, i))
enddo
Omega_x=0d0;Omega_y=0d0; Omega_z=0d0
!do m= NumOccupied , NumOccupied
! do n= NumOccupied+1, NumOccupied+1
do m= 1, NumOccupied
do n= 1, Num_wann
if (m==n) cycle
Omega_x(m)= Omega_x(m)+ velocity_Ham(n, m, 2)*velocity_Ham(m, n, 3)/((W(m)-W(n))**2)
Omega_y(m)= Omega_y(m)+ velocity_Ham(n, m, 3)*velocity_Ham(m, n, 1)/((W(m)-W(n))**2)
Omega_z(m)= Omega_z(m)+ velocity_Ham(n, m, 1)*velocity_Ham(m, n, 2)/((W(m)-W(n))**2)
enddo ! m
enddo ! n
Omega_x= -Omega_x*2d0*zi
Omega_y= -Omega_y*2d0*zi
Omega_z= -Omega_z*2d0*zi
return
end subroutine berry_curvarture_singlek_numoccupied
subroutine berry_curvarture_singlek_allbands(Dmn_Ham, Omega_BerryCurv)
!> Calculate Berry curvature for a sigle k point and all bands
!> ref : eqn (30) Physical Review B 74, 195118(2006)
!> \Omega_n^{\gamma}(k)=i\sum_{\alpha\beta}\epsilon_{\gamma\alpha\beta}(D^{\alpha\dag}D^{\beta})_{nn}
!> Dec. 30 2019 by Quansheng Wu
!> Copyright (c) 2019 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
!> D_mn^H=V_mn/(En-Em) for m!=n
!> D_nn^H=0
complex(dp), intent(in) :: Dmn_Ham(Num_wann, Num_wann, 3)
!> Berry curvature vectors for all bands
real(dp), intent(out) :: Omega_BerryCurv(Num_wann, 3)
integer :: m, n, i, alphai(3), betai(3)
complex(dp), allocatable :: Amat(:, :)
complex(dp), allocatable :: Dmn_Ham_dag(:, :, :)
complex(dp) :: zz(3)
allocate(Amat(Num_wann, Num_wann))
allocate(Dmn_Ham_dag(Num_wann, Num_wann, 3))
Amat= 0d0
do i=1, 3
Dmn_Ham_dag(:, :, i)=conjg(transpose(Dmn_Ham(:, :, i)))
! Dmn_Ham_dag(:, :, i)=-(Dmn_Ham(:, :, i))
enddo
!print *, Dmn_Ham_dag(1, 2, 2)
!print *, Dmn_Ham(2, 1, 3)
Omega_BerryCurv= 0d0
alphai(1)=2;alphai(2)=3;alphai(3)=1;
betai(1) =3; betai(2)=1; betai(3)=2;
zz=0d0
do i=1, 3
call mat_mul(Num_wann, Dmn_Ham_dag(:, :, alphai(i)), Dmn_Ham(:, :, betai(i)), Amat)
do m=1, Num_wann
!Omega_BerryCurv(m, i)= Omega_BerryCurv(m, i)- aimag(Amat(m, m))
do n=1, Num_wann
Omega_BerryCurv(m, i)= Omega_BerryCurv(m, i)- &
aimag(Dmn_Ham_dag(m, n, alphai(i))*Dmn_Ham(n, m, betai(i)))
enddo
enddo
call mat_mul(Num_wann, Dmn_Ham_dag(:, :, betai(i)), Dmn_Ham(:, :, alphai(i)), Amat)
do m=1, Num_wann
!Omega_BerryCurv(m, i)= Omega_BerryCurv(m, i)+ aimag(Amat(m, m))
do n=1, Num_wann
Omega_BerryCurv(m, i)= Omega_BerryCurv(m, i)+ &
aimag(Dmn_Ham_dag(m, n, betai(i))*Dmn_Ham(n, m, alphai(i)))
enddo
enddo
enddo
return
end subroutine berry_curvarture_singlek_allbands
subroutine orbital_magnetization_singlek_allbands(Dmn_Ham, Vmn_Ham_nondiag, m_OrbMag)
!> Calculate orbital magnetization for a sigle k point and all bands
!> m_n^{\gamma}(k)=-i*e/2/\hbar\sum_{\alpha\beta}\epsilon_{\gamma\alpha\beta}(D^{\alpha\dag}V^{\beta})_{nn}
!> Dec. 30 2019 by Quansheng Wu
!> Copyright (c) 2019 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
!> D_mn^H=V_mn/(En-Em) for m!=n
!> D_nn^H=0
complex(dp), intent(in) :: Dmn_Ham(Num_wann, Num_wann, 3)
!> Vmn_Ham_nondiag(m, n, i)= Vmn_Ham(m, n, i) for m!=n
!> Vmn_Ham_nondiag(n, n, i)= 0
complex(dp), intent(in) :: Vmn_Ham_nondiag(Num_wann, Num_wann, 3)
!> Berry curvature vectors for all bands
real(dp), intent(out) :: m_OrbMag(Num_wann, 3)
integer :: m, i, ialpha(3), ibeta(3)
complex(dp), allocatable :: Amat(:, :)
complex(dp), allocatable :: Dmn_Ham_dag(:, :, :)
allocate(Amat(Num_wann, Num_wann))
allocate(Dmn_Ham_dag(Num_wann, Num_wann, 3))
Amat=0d0
Dmn_Ham_dag = -Dmn_Ham
m_OrbMag= 0d0
ialpha(1)=2;ialpha(2)=3;ialpha(3)=1;
ibeta(1)=3;ibeta(2)=1;ibeta(3)=2;
do i=1, 3
call mat_mul(Num_wann, Dmn_Ham_dag(:, :, ialpha(i)), Vmn_Ham_nondiag(:, :, ibeta(i)), Amat)
do m=1, Num_wann
m_OrbMag(m, i)= m_OrbMag(m, i)+ aimag(Amat(m, m))
enddo
call mat_mul(Num_wann, Dmn_Ham_dag(:, :, ibeta(i)), Vmn_Ham_nondiag(:, :, ialpha(i)), Amat)
do m=1, Num_wann
m_OrbMag(m, i)= m_OrbMag(m, i)- aimag(Amat(m, m))
enddo
enddo
!> 1/2d0 comes from e/\hbar/2d0, here we use atomic unit
m_OrbMag= m_OrbMag/2d0
!> if we want to use Bohr magneton \mu_B as unit, we need multiply it by two.
return
end subroutine orbital_magnetization_singlek_allbands
subroutine get_Vmn_Ham_nondiag(velocity_Ham, velocity_Ham_nondiag)
!> the diagonal of velocity_Ham_nondiag is zero
!> and velocity_Ham_nondiag(m, n)= -velocity_Ham(m, n)
use para, only : dp, Num_wann
implicit none
complex(dp), intent(in) :: velocity_Ham(Num_wann, Num_wann, 3)
complex(dp), intent(out) :: velocity_Ham_nondiag(Num_wann, Num_wann, 3)
integer :: n
velocity_Ham_nondiag= -velocity_Ham
do n=1, Num_wann
velocity_Ham_nondiag(n, n, :)=0d0
enddo
return
end subroutine get_Vmn_Ham_nondiag
subroutine get_Dmn_Ham(W, velocity_Ham, Dmn_Ham)
!> define D matrix
!> Eq (24) Phys.Rev.B 74, 195118(2006)
!> Dmn_Ham_dag=-Dmn_Ham
use para, only : dp, Num_wann, eps12, zzero
implicit none
real(dp), intent(in) :: W(Num_wann)
complex(dp), intent(in) :: velocity_Ham(Num_wann, Num_wann, 3)
complex(dp), intent(out) :: Dmn_Ham(Num_wann, Num_wann, 3)
integer :: m, n, i
Dmn_Ham = zzero
do i=1, 3
do n=1, Num_wann
do m=1, Num_wann
if (abs(W(n)-W(m))>eps12) then
Dmn_Ham(m,n,i)=velocity_Ham(m,n,i)/(W(n)-W(m))
endif
enddo
Dmn_Ham(n, n, i)=zzero
enddo
enddo
return
end subroutine get_Dmn_Ham
subroutine berry_curvarture_slab
!> Calculate Berry curvature
!
!> ref : Physical Review B 74, 195118(2006)
!
!> Aug. 06 2018 by Quansheng Wu @ EPFL
!
! Copyright (c) 2018 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
integer :: ik, ierr, Mdim, i, j
real(dp) :: k(2)
!> k points slice
real(dp), allocatable :: k12(:, :)
real(dp), allocatable :: k12_xyz(:, :)
real(dp), external :: norm
!> Berry curvature (k)
complex(dp), allocatable :: Omega_z(:)
complex(dp), allocatable :: Omega(:)
complex(dp), allocatable :: Omega_mpi(:)
Mdim = Num_wann*Nslab
allocate( k12(2, Nk1*Nk2))
allocate( k12_xyz(2, Nk1*Nk2))
allocate( Omega_z(Mdim))
allocate( Omega (Nk1*Nk2))
allocate( Omega_mpi(Nk1*Nk2))
k12=0d0
k12_xyz=0d0
omega= 0d0
omega_mpi= 0d0
!> k12 is centered at K2d_start
ik=0
do i= 1, nk1
do j= 1, nk2
ik=ik+1
k12(:, ik)=K2D_start+ (i-1)*K2D_vec1/dble(nk1-1) &
+ (j-1)*K2D_vec2/dble(nk2-1)- (K2D_vec1+K2D_vec2)/2d0
k12_xyz(:, ik)= k12(1, ik)* Ka2+ k12(2, ik)* Kb2
enddo
enddo
do ik= 1+ cpuid, Nk1*Nk2, num_cpu
if (cpuid==0) write(stdout, *)'Berry curvature ik, nk1*nk2 ', ik, Nk1*Nk2
!> diagonalize hamiltonian
k= k12(:, ik)
Omega_z= 0d0
call berry_curvarture_singlek_numoccupied_slab_total(k, Omega_z(1))
Omega(ik) = sum(Omega_z)
enddo ! ik
Omega_mpi= 0d0
#if defined (MPI)
call mpi_allreduce(Omega,Omega_mpi,size(Omega_mpi),&
mpi_dc,mpi_sum,mpi_cmw,ierr)
#else
Omega_mpi= Omega
#endif
!> output the Berry curvature to file
outfileindex= outfileindex+ 1
if (cpuid==0) then
open(unit=outfileindex, file='Berrycurvature_slab.dat')
write(outfileindex, '(20a28)')'# Unit of Berry curvature is Angstrom^2'
write(outfileindex, '(20a28)')'# kx (1/A)', 'ky (1/A)', &
'Omega_z'
ik= 0
do i= 1, nk1
do j= 1, nk2
ik= ik+ 1
write(outfileindex, '(20E28.10)')k12_xyz(:, ik)*Angstrom2atomic, real(Omega_mpi(ik))/Angstrom2atomic**2
enddo
write(outfileindex, *) ' '
enddo
close(outfileindex)
endif
!> generate gnuplot script to plot the Berry curvature
outfileindex= outfileindex+ 1
if (cpuid==0) then
open(unit=outfileindex, file='Berrycurvature_slab.gnu')
write(outfileindex, '(a)')"set encoding iso_8859_1"
write(outfileindex, '(a)')'#set terminal pngcairo truecolor enhanced size 1920, 1680 font ",60"'
write(outfileindex, '(a)')'set terminal png truecolor enhanced size 1920, 1680 font ",60"'
write(outfileindex, '(a)')"set output 'Berrycurvature_slab.png'"
write(outfileindex, '(a)')"set palette rgbformulae 33,13,10"
write(outfileindex, '(a)')"unset ztics"
write(outfileindex, '(a)')"set size ratio -1"
write(outfileindex, '(a)')"set size 0.9, 0.95"
write(outfileindex, '(a)')"set origin 0.05, 0.02"
write(outfileindex, '(a)')"unset key"
write(outfileindex, '(a)')"set pm3d"
write(outfileindex, '(a)')"#set zbrange [ -10: 10] "
write(outfileindex, '(a)')"#set cbrange [ -100: 100] "
write(outfileindex, '(a)')"set view map"
write(outfileindex, '(a)')"set border lw 3"
write(outfileindex, '(a)')"set xlabel 'k (1/{\305})'"
write(outfileindex, '(a)')"set ylabel 'k (1/{\305})'"
write(outfileindex, '(a)')"set xrange [] noextend"
write(outfileindex, '(a)')"set yrange [] noextend"
write(outfileindex, '(a)')"set ytics 0.5 nomirror scale 0.5"
write(outfileindex, '(a)')"set pm3d interpolate 2,2"
write(outfileindex, '(a)')"set title 'Berry Curvature {/Symbol W}_z (Ang^2)'"
write(outfileindex, '(a)')"set colorbox"
write(outfileindex, '(a)')"splot 'Berrycurvature_slab.dat' u 1:2:3 w pm3d"
close(outfileindex)
endif
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
deallocate( k12)
deallocate( k12_xyz)
deallocate( Omega_z)
deallocate( Omega, Omega_mpi)
return
end subroutine berry_curvarture_slab
subroutine berry_curvarture_line
!> Calculate Berry curvature for a k line defined bu KPATH_BULK
!
!> ref : Physical Review B 74, 195118(2006)
!
!> eqn (9), Eqn (34)
!
!> Sep. 28 2018 by Quansheng Wu @ EPFL
!
! Copyright (c) 2018 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
integer :: ik, ierr, i
real(dp) :: k(3), o1(3), k_cart(3)
real(dp) :: time_start, time_end, time_start0, ybound_min, ybound_max
real(dp), external :: norm
!> Berry curvature (3, bands, k)
complex(dp), allocatable :: Omega_x(:), Omega_y(:), Omega_z(:)
real(dp), allocatable :: Omega(:, :), Omega_mpi(:, :)
!> energy bands
real(dp), allocatable :: eigv(:,:)
real(dp), allocatable :: eigv_mpi(:,:)
allocate( Omega_x(Num_wann))
allocate( Omega_y(Num_wann))
allocate( Omega_z(Num_wann))
allocate( eigv (Num_wann, nk3_band))
allocate( eigv_mpi(Num_wann, nk3_band))
allocate( Omega (3, nk3_band))
allocate( Omega_mpi(3, nk3_band))
omega= 0d0
omega_mpi= 0d0
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
if (cpuid==0.and. mod(ik/num_cpu, 100)==0) &
write(stdout, '(a, i9, " /", i10, a, f10.1, "s", a, f10.1, "s")') &
' Berry curvature: ik', ik, nk3_band, ' time left', &
(nk3_band-ik)*(time_end- time_start)/num_cpu, &
' time elapsed: ', time_end-time_start0
!> a k point in fractional coordinates
k= k3points(:, ik)
call now(time_start)
Omega_x= 0d0
Omega_y= 0d0
Omega_z= 0d0
!call berry_curvarture_singlek_numoccupied_old(k, Omega_x, Omega_y, Omega_z)
if (Berrycurvature_kpath_EF_calc) then
call berry_curvarture_singlek_EF(k, E_arc, Omega_x, Omega_y, Omega_z)
else if (BerryCurvature_kpath_Occupied_calc) then
call berry_curvarture_singlek_numoccupied_total(k, Omega_x(1), Omega_y(1), Omega_z(1))
else
write(*, *) 'ERROR: In subroutine berry_curvarture_line, we only support BerryCurvature_kpath_Occupied_calc '
write(*, *) ' and Berrycurvature_kpath_EF_calc '
stop
endif
Omega(1, ik) = real(sum(Omega_x))
Omega(2, ik) = real(sum(Omega_y))
Omega(3, ik) = real(sum(Omega_z))
call now(time_end)
enddo ! ik
Omega_mpi= 0d0
#if defined (MPI)
call mpi_allreduce(Omega,Omega_mpi,size(Omega_mpi),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
#else
Omega_mpi= Omega
#endif
!> output the Berry curvature to file
outfileindex= outfileindex+ 1
if (cpuid==0) then
open(unit=outfileindex, file='Berrycurvature_line.dat')
write(outfileindex, '(20a18)')'# Column 1 kpath with accumulated length in the kpath'
write(outfileindex, '(20a18)')'# Column 2-4 Berry curvature (Ang^2)'
write(outfileindex, '(20a18)')'# k (1/A)', &
'real(Omega_x)', 'real(Omega_y)', 'real(Omega_z)'
do ik= 1, nk3_band
k=k3points(:, ik)
write(outfileindex, '(20E18.8)')k3len(ik)*Angstrom2atomic, real(Omega_mpi(:, ik))/Angstrom2atomic**2
enddo
close(outfileindex)
endif
ybound_min= minval(real(Omega_mpi))-2
ybound_max= maxval(real(Omega_mpi))+5
!> write script for gnuplot
outfileindex= outfileindex+ 1
if (cpuid==0) then
open(unit=outfileindex, file='Berrycurvature_line.gnu')
write(outfileindex, '(a)')"set encoding iso_8859_1"
write(outfileindex, '(a)') 'set terminal pdf enhanced color font ",16"'
write(outfileindex, '(a)')"set output 'Berrycurvature_line.pdf'"
write(outfileindex, '(a)')'set style data linespoints'
write(outfileindex, '(a)')'unset key'
write(outfileindex, '(a)')'set pointsize 0.8'
write(outfileindex, '(a)')'set view 0,0'
write(outfileindex, '(a, f10.5, a)')'set xrange [0: ', maxval(k3len)*Angstrom2atomic, ']'
if (index(Particle,'phonon')/=0) then
write(outfileindex, '(a, f10.5, a)')'set yrange [0:', ybound_max, ']'
write(outfileindex, '(a)')'set ylabel "Frequency (THz)"'
else
write(outfileindex, '(a)')'set ylabel "Energy (eV)"'
write(outfileindex, '(a, f10.5, a, f10.5, a)')'set yrange [', ybound_min, ':', ybound_max, ']'
endif
write(outfileindex, 202, advance="no") (k3line_name(i), k3line_stop(i)*Angstrom2atomic, i=1, nk3lines)
write(outfileindex, 203)k3line_name(nk3lines+1), k3line_stop(nk3lines+1)*Angstrom2atomic
do i=1, nk3lines-1
if (index(Particle,'phonon')/=0) then
write(outfileindex, 204)k3line_stop(i+1)*Angstrom2atomic, 0.0, k3line_stop(i+1)*Angstrom2atomic, ybound_max
else
write(outfileindex, 204)k3line_stop(i+1)*Angstrom2atomic, ybound_min, k3line_stop(i+1)*Angstrom2atomic, ybound_max
endif
enddo
write(outfileindex, '(a)')"plot 'Berrycurvature_line.dat' \"
write(outfileindex, '(a)')"u 1:2 w lp lc rgb 'red' lw 2 pt 7 ps 0.2 title '{/Symbol W}_x', \"
write(outfileindex, '(a)')"'' u 1:3 w lp lc rgb 'green' lw 2 pt 7 ps 0.2 title '{/Symbol W}_y', \"
write(outfileindex, '(a)')"'' u 1:4 w lp lc rgb 'blue' lw 2 pt 7 ps 0.2 title '{/Symbol W}_z' "
close(outfileindex)
endif
202 format('set xtics (',20('"',A3,'" ',F10.5,','))
203 format(A3,'" ',F10.5,')')
204 format('set arrow from ',F10.5,',',F10.5,' to ',F10.5,',',F10.5, ' nohead')
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
deallocate( Omega_x, Omega_y, Omega_z)
deallocate( Omega, Omega_mpi)
return
end subroutine berry_curvarture_line
subroutine berry_curvarture_cube
!> Calculate Berry curvature in a cube defined in KCUBE_BULK
!
!> ref : Physical Review B 74, 195118(2006)
!
!> eqn (34)
!
!> Sep. 28 2018 by Quansheng Wu @ EPFL
!
! Copyright (c) 2018 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
integer :: ik, ierr, ikx, iky, ikz, knv3, i, m, n
real(dp) :: k(3), o1(3), k_cart(3)
real(dp) :: time_start, time_end, time_start0
real(dp), external :: norm
!> Berry curvature and orbital magnetization (3, bands, k)
real(dp), allocatable :: Omega_allk(:, :, :), Omega_allk_mpi(:, :, :)
real(dp), allocatable :: m_OrbMag_allk(:, :, :), m_OrbMag_allk_mpi(:, :, :)
real(dp), allocatable :: Omega_BerryCurv(:, :), m_OrbMag(:, :)
!> velocity matrix in Wannier basis
complex(dp), allocatable :: Vmn_wann(:, :, :)
!> velocity matrix in Hamiltonian basis
complex(dp), allocatable :: Vmn_Ham(:, :, :), Vmn_Ham_nondiag(:, :, :)
complex(dp), allocatable :: Dmn_Ham(:, :, :)
!> Hamiltonian, eigenvalue and eigenvectors
complex(dp), allocatable :: UU(:, :)
real(dp), allocatable :: W(:)
real(dp), allocatable :: eigval_allk(:, :), eigval_allk_mpi(:, :)
knv3=Nk1*Nk2*Nk3
allocate(Vmn_wann(Num_wann, Num_wann, 3), Vmn_Ham(Num_wann, Num_wann, 3))
allocate(Dmn_Ham(Num_wann,Num_wann,3), Vmn_Ham_nondiag(Num_wann, Num_wann, 3))
allocate(W(Num_wann))
allocate(UU(Num_wann, Num_wann))
allocate(eigval_allk(Num_wann, knv3))
allocate(eigval_allk_mpi(Num_wann, knv3))
allocate( Omega_BerryCurv(Num_wann, 3), m_OrbMag(Num_wann, 3))
allocate( Omega_allk (Num_wann, 3, knv3))
allocate( Omega_allk_mpi(Num_wann, 3, knv3))
allocate( m_OrbMag_allk (Num_wann, 3, knv3))
allocate( m_OrbMag_allk_mpi(Num_wann, 3, knv3))
Omega_BerryCurv= 0d0
m_OrbMag=0d0
Omega_allk= 0d0
eigval_allk= 0d0
m_OrbMag_allk=0d0
Omega_allk_mpi= 0d0
eigval_allk_mpi= 0d0
m_OrbMag_allk_mpi=0d0
time_start= 0d0
time_start0= 0d0
call now(time_start0)
time_start= time_start0
time_end = time_start0
do ik= 1+ cpuid, knv3, num_cpu
if (cpuid==0.and. mod(ik/num_cpu, 100)==0) &
write(stdout, '(a, i9, " /", i10, a, f10.1, "s", a, f10.1, "s")') &
' Berry curvature: ik', ik, knv3, ' time left', &
(knv3-ik)*(time_end- time_start)/num_cpu, &
' time elapsed: ', time_end-time_start0
ikx= (ik-1)/(nk2*nk3)+1
iky= ((ik-1-(ikx-1)*Nk2*Nk3)/nk3)+1
ikz= (ik-(iky-1)*Nk3- (ikx-1)*Nk2*Nk3)
k= K3D_start_cube+ K3D_vec1_cube*(ikx-1)/dble(nk1) &
+ K3D_vec2_cube*(iky-1)/dble(nk2) &
+ K3D_vec3_cube*(ikz-1)/dble(nk3)
!- (K3D_vec1_cube+ K3D_vec2_cube+ K3D_vec3_cube)/2d0
call now(time_start)
!> calculation bulk hamiltonian by a direct Fourier transformation of HmnR
call ham_bulk_atomicgauge(k, UU)
!> diagonalization by call zheev in lapack
call eigensystem_c( 'V', 'U', Num_wann, UU, W)
eigval_allk(:, ik) = W
!> get velocity operator in Hamiltonian basis
call dHdk_atomicgauge_Ham(k, UU, Vmn_Ham)
call get_Dmn_Ham(W, Vmn_Ham, Dmn_Ham)
call get_Vmn_Ham_nondiag(Vmn_Ham, Vmn_Ham_nondiag)
call berry_curvarture_singlek_allbands(Dmn_Ham, Omega_BerryCurv)
call orbital_magnetization_singlek_allbands(Dmn_Ham, Vmn_Ham_nondiag, m_OrbMag)
Omega_allk(:, :, ik) = Omega_BerryCurv
m_OrbMag_allk(:, :, ik) = m_OrbMag
call now(time_end)
enddo ! ik
#if defined (MPI)
call mpi_allreduce(Omega_allk,Omega_allk_mpi,size(Omega_allk_mpi),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
call mpi_allreduce(eigval_allk,eigval_allk_mpi,size(eigval_allk_mpi),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
call mpi_allreduce(m_OrbMag_allk,m_OrbMag_allk_mpi,size(m_OrbMag_allk_mpi),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
#else
Omega_allk_mpi= Omega_allk
eigval_allk_mpi= eigval_allk
m_OrbMag_allk_mpi= m_OrbMag_allk
#endif
!> write out Berry curvature and orbital magnetization to a file which
!> can be open by software Fermisurfer.
outfileindex= outfileindex+ 1
if (cpuid==0) then
open(unit=outfileindex, file='Berrycurvature_norm.frmsf')
write(outfileindex, *) Nk1, Nk2, Nk3
write(outfileindex, *) 1
write(outfileindex, *) Num_wann
write(outfileindex, '(3f12.6)') Origin_cell%Kua
write(outfileindex, '(3f12.6)') Origin_cell%Kub
write(outfileindex, '(3f12.6)') Origin_cell%Kuc
do m=1, Num_wann
do ik= 1, knv3
write(outfileindex, '(E18.10)') eigval_allk_mpi(m, ik)-E_arc
enddo
enddo
do m=1, Num_wann
do ik= 1, knv3
o1= Omega_allk_mpi(m, :, ik)/Angstrom2atomic**2
write(outfileindex, '(E18.10)') norm(o1)
enddo
enddo
close(outfileindex)
endif
!> write out Berry curvature and orbital magnetization to a file which
!> can be open by software Fermisurfer.
outfileindex= outfileindex+ 1
if (cpuid==0) then
open(unit=outfileindex, file='Berrycurvature_x.frmsf')
write(outfileindex, *) Nk1, Nk2, Nk3
write(outfileindex, *) 1
write(outfileindex, *) Num_wann
write(outfileindex, '(3f12.6)') Origin_cell%Kua
write(outfileindex, '(3f12.6)') Origin_cell%Kub
write(outfileindex, '(3f12.6)') Origin_cell%Kuc
do m=1, Num_wann
do ik= 1, knv3
write(outfileindex, '(E18.10)') eigval_allk_mpi(m, ik)-E_arc
enddo
enddo
do m=1, Num_wann
do ik= 1, knv3
o1= Omega_allk_mpi(m, :, ik)/Angstrom2atomic**2
write(outfileindex, '(E18.10)') o1(1)
enddo
enddo
close(outfileindex)
endif
!> write out Berry curvature and orbital magnetization to a file which
!> can be open by software Fermisurfer.
outfileindex= outfileindex+ 1
if (cpuid==0) then
open(unit=outfileindex, file='Berrycurvature_z.frmsf')
write(outfileindex, *) Nk1, Nk2, Nk3
write(outfileindex, *) 1
write(outfileindex, *) Num_wann
write(outfileindex, '(3f12.6)') Origin_cell%Kua
write(outfileindex, '(3f12.6)') Origin_cell%Kub
write(outfileindex, '(3f12.6)') Origin_cell%Kuc
do m=1, Num_wann
do ik= 1, knv3
write(outfileindex, '(E18.10)') eigval_allk_mpi(m, ik)-E_arc
enddo
enddo
do m=1, Num_wann
do ik= 1, knv3
o1= Omega_allk_mpi(m, :, ik)/Angstrom2atomic**2
write(outfileindex, '(E18.10)') o1(3)
enddo
enddo
close(outfileindex)
endif
!> orbital_magnetization
outfileindex= outfileindex+ 1
if (cpuid==0) then
open(unit=outfileindex, file='Orbital_magnetization_norm.frmsf')
write(outfileindex, *) Nk1, Nk2, Nk3
write(outfileindex, *) 1
write(outfileindex, *) Num_wann
write(outfileindex, '(3f12.6)') Origin_cell%Kua
write(outfileindex, '(3f12.6)') Origin_cell%Kub
write(outfileindex, '(3f12.6)') Origin_cell%Kuc
do m=1, Num_wann
do ik= 1, knv3
write(outfileindex, '(E18.10)') eigval_allk_mpi(m, ik)-E_arc
enddo
enddo
do m=1, Num_wann
do ik= 1, knv3
o1= m_OrbMag_allk_mpi(m, :, ik)
write(outfileindex, '(E18.10)') norm(o1)
enddo
enddo
close(outfileindex)
endif
!> orbital_magnetization
outfileindex= outfileindex+ 1
if (cpuid==0) then
open(unit=outfileindex, file='Orbital_magnetization_z.frmsf')
write(outfileindex, *) Nk1, Nk2, Nk3
write(outfileindex, *) 1
write(outfileindex, *) Num_wann
write(outfileindex, '(3f12.6)') Origin_cell%Kua
write(outfileindex, '(3f12.6)') Origin_cell%Kub
write(outfileindex, '(3f12.6)') Origin_cell%Kuc
do m=1, Num_wann
do ik= 1, knv3
write(outfileindex, '(E18.10)') eigval_allk_mpi(m, ik)-E_arc
enddo
enddo
do m=1, Num_wann
do ik= 1, knv3
o1= m_OrbMag_allk_mpi(m, :, ik)
write(outfileindex, '(E18.10)') o1(3)