-
Notifications
You must be signed in to change notification settings - Fork 0
/
sigma_OHE.f90
executable file
·3199 lines (2759 loc) · 120 KB
/
sigma_OHE.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
!---------------------------------------------------------!
!> Calculate magnetoresistance with R.G.Chambers's formula!
!> based on Boltzmann transport !
!> Written By QuanSheng Wu ([email protected]) !
!> Modified by Hanqi Pi ([email protected]) !
!> Thanks Yi Liu for helpful discussions !
!> !
!> References : !
!> [1] Electrons in metals and semiconductors, !
!> R.G. Chambers, !
!> [2] Ab initio investigation of magnetic transport !
!> properties by Wannier interpolation, PHYSICAL !
!> REVIEW B 79, 245123 (2009), Yi Liu, Hai-Jun Zhang, !
!> and Yugui Yao !
!> [3] Magnetoresistance from Fermi surface topology, !
!> ShengNan Zhang, QuanSheng Wu, Yi Liu, and !
!> Oleg V. Yazyev,Phys. Rev. B 99, 035142 (2019) !
!> !
!> Implemented on Oct. 07 2017 !
!> uploaded on Sep. 05. 2019 !
!---------------------------------------------------------!
subroutine sigma_ohe_calc_symm(mu_array, KBT_array, BTau_array, Nband_Fermi_Level, bands_fermi_level, sigma_ohe_tensor)
!-----------------------------------------------------------!
!> This is the main routine to calculate the !
!> magnetoresistance/conductivity based on Boltzmann !
!> transport equation. !
!> !
!> We calculate the conductivity and resistivity under !
!> the band-resolved constantly relaxation time. We give !
!> the simga/tau of Btau and rho*tau of Btau instead of !
!> sigma and rho. !
!> !
!> This version impose symmetry constrain to get the !
!> transport coefficients under magnetic field. !
!> !
!> It produces 4 types of output files, !
!> every type has OmegaNum files at different chemical !
!> !
!> 1. sigma_bands_mu_***eV.dat : the band-resolved sigma/tau!
!> 2. rho_bands_mu_***eV.dat : the band-resolved rho*tau !
!> 3. sigma_total_mu_***eV.dat : the total sigma/tau !
!> 4. rho_total_mu_***eV.dat : the total rho*tau !
!> !
!> The first two files are output on a grid of !
!> (Nband_Fermi_Level, NumT, NBTau), while the last two !
!> files are output on a grid of (NumT, NBTau). NumT is the !
!> number of temperature points, NBTau is the number of Btau!
!> and Nband_Fermi_Level is the number of bands. !
!> !
!> The output sigma/tau is in the unit of 1/(Ohm*m*s), !
!> while the output rho*tau is in the unit of Ohm*m*s. !
!-----------------------------------------------------------!
use wmpi
use para
implicit none
integer, intent(inout) :: Nband_Fermi_Level
integer, intent(inout) :: bands_fermi_level(Nband_Fermi_Level)
real(dp), intent(in) :: KBT_array(NumT) ! K_Boltzmann*Temperature in eV
real(dp), intent(in) :: mu_array(OmegaNum) ! chemical potential relative to Fermi level in eV
real(dp), intent(in) :: BTau_array(NBTau) ! omega*tau without units
real(dp), intent(inout) :: sigma_ohe_tensor(9, NBTau, OmegaNum, NumT, Nband_Fermi_Level)
real(dp) :: coeff, mu, BTau, KBT, exponent_max
integer :: ie, ibtau, ikt
integer :: Nk_total, Nk_current, Nk_start, Nk_end
integer :: knv3_left, knv3_left_mod, knv3, knv3_mod
integer :: ik, iband, ik1, ik2, ik3, ik_temp
integer :: ierr, it, i, ix, j1, j2, j
integer :: nrecevs
real(dp) :: v_t(3), vo_k, v_k(3)
real(dp) :: k(3), k_start(3), magnetic_field(3)
real(dp) :: sigma_symm_t(9), rho(3, 3)
real(dp) :: time_start, time_end
real(dp) :: time_start0, time_end0
integer :: NSlice_Btau_inuse
!> Btau slices for Runge-Kutta integration
real(dp) :: Btau_start, Btau_final, DeltaBtau
logical :: fail
integer :: icycle, ishift
!> energy bands
real(dp) :: EE
real(dp), allocatable :: Ek(:) ! in eV
real(dp), allocatable :: Enk(:, :)
!> minus fermi derivation
real(dp) :: minusdfde
!> 3-component velocity for each band and each k point
real(dp), allocatable :: velocity_k(:, :)
real(dp), allocatable :: velocity_bar_k(:)
real(dp), allocatable :: velocity(:, :, :)
!> fileindex list for sigma_bands and rho_bands
integer, allocatable :: sigmafileindex(:)
integer, allocatable :: rhofileindex(:)
!> 3-component velocity for each band and each k point
!> Eq.(3) in PRB 79, 245123(2009)
real(dp), allocatable :: velocity_bar(:, :, :)
type(kcube_type) :: KCube3D_total
type(kcube_type) :: KCube3D_left(Nband_Fermi_Level)
!> some arrays for mpi
integer, allocatable :: info(:, :) !> howmany integers to be sent on each cpu
integer, allocatable :: Displs(:, :)
!> number of steps used in the Runge-Kutta integration
integer :: NSlice_Btau
integer :: NSlice_Btau_local
real(dp), allocatable :: kout(:, :)
!> define some arrays for different bands. Since there are different number
!> of k points left for different bands.
type klist_iband_type
!> dim=3*NSlice_Btau
real(dp), allocatable :: klist_rkfs(:, :)
!> dim=3*NSlice_Btau
real(dp), allocatable :: velocity_k(:, :)
!> calculate -df(e)/de, where f(e) is the Fermi-Dirac distribution
!> dim=(NumT, OmegaNum))
real(dp), allocatable :: minusdfde(:, :)
end type klist_iband_type
type(klist_iband_type) :: klist_iband(Nband_Fermi_Level)
type sigma_iband_type
! sigma_ohe_tensor_k(9, NBTau, OmegaNum, NumT)
real(dp), allocatable :: sigma_ohe_tensor_k_mpi(:, :, :, :)
real(dp), allocatable :: sigma_ohe_tensor_k(:, :, :, :)
! sigma_ohe_tensor_kz(9, NBTau, OmegaNum, NumT, Nk3)
! real(dp), allocatable :: sigma_ohe_tensor_kz(:, :, :, :, :)
!> time cost for k point
real(dp), allocatable :: time_cost(:)
real(dp), allocatable :: time_cost_mpi(:)
end type sigma_iband_type
type(sigma_iband_type) :: sigma_iband_k(Nband_Fermi_Level)
!> file index
!integer, allocatable :: myfileindex(:)
character(80) :: sigmafilename, bandname, tname, muname
!> inverse of group operator
real(dp) :: Tmat(3, 3)
real(dp), allocatable :: pgop_cart_inverse(:, :, :)
real(dp), external :: det3
! default value
icycle = 1
fail = .False.
allocate(pgop_cart_inverse(3, 3, number_group_operators))
do i= 1, number_group_operators
Tmat= pgop_cart(:, :, i)
call inv_r(3, Tmat)
pgop_cart_inverse(:, :, i)= Tmat
enddo
!> set up
!> Distribute k kpoints in the 3DCube into different MPI threads
knv3= KCube3D_symm%Nk_total_symm
KCube3D_total%Nk_total= knv3
knv3_mod= mod(knv3, num_cpu)
if (knv3_mod==0) then !> perfect divided
KCube3D_total%Nk_current= knv3/num_cpu
KCube3D_total%Nk_start=1+ knv3*cpuid/num_cpu
KCube3D_total%Nk_end =(1+cpuid)*knv3/num_cpu
else if (knv3/num_cpu==0) then !> Number of MPI threads is large than knv3
KCube3D_total%Nk_current= 1 !> one k piont per MPI thread
KCube3D_total%Nk_start= cpuid+ 1 !> one k piont per MPI thread
KCube3D_total%Nk_end = cpuid+ 1
if (cpuid+1 > knv3) then
KCube3D_total%Nk_start= 1
KCube3D_total%Nk_end = 0
endif
else
KCube3D_total%Nk_current= knv3/num_cpu+ 1
if (cpuid< knv3_mod) then
KCube3D_total%Nk_start= 1+ cpuid*KCube3D_total%Nk_current
KCube3D_total%Nk_end = (1+cpuid)*KCube3D_total%Nk_current
else
KCube3D_total%Nk_start= knv3_mod*KCube3D_total%Nk_current+ &
(cpuid-knv3_mod)*(KCube3D_total%Nk_current-1)+1
KCube3D_total%Nk_end = knv3_mod*KCube3D_total%Nk_current+ &
(cpuid-knv3_mod+1)*(KCube3D_total%Nk_current-1)
endif
endif
!> calculate the volume of the k cube
allocate(KCube3D_total%k_direct(3, KCube3D_total%Nk_start:KCube3D_total%Nk_end))
allocate(KCube3D_total%weight_k(KCube3D_total%Nk_start:KCube3D_total%Nk_end))
do ik= KCube3D_total%Nk_start, KCube3D_total%Nk_end
j1= KCube3D_symm%ik_array_symm(ik)
ik1= (j1-1)/(Nk2*Nk3)+1
ik2= ((j1-1-(ik1-1)*Nk2*Nk3)/Nk3)+1
ik3= (j1-(ik2-1)*Nk3- (ik1-1)*Nk2*Nk3)
KCube3D_total%k_direct(1, ik)= (ik1-1d0)/dble(Nk1)
KCube3D_total%k_direct(2, ik)= (ik2-1d0)/dble(Nk2)
KCube3D_total%k_direct(3, ik)= (ik3-1d0)/dble(Nk3)
KCube3D_total%weight_k(ik)= KCube3D_symm%weight_k(ik)
enddo
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
!> setup NSlice_Btau
!> NSlice_Btau should be the integer times of NBTau
if (NBTau>1) then
NSlice_Btau= (NBTau-1)*(Nslice_BTau_Max/(NBTau-1))
else
NSlice_Btau= 1
endif
if (cpuid.eq.0) write(stdout, *) ' NSlice_Btau :', NSlice_Btau
Nk_total= KCube3D_total%Nk_total
Nk_current= KCube3D_total%Nk_current
Nk_start= KCube3D_total%Nk_start
Nk_end= KCube3D_total%Nk_end
allocate( Ek(Nband_Fermi_Level))
allocate( Enk(Nk_start:Nk_end, Nband_Fermi_Level))
allocate( velocity(3, Nk_start:Nk_end, Nband_Fermi_Level))
allocate( velocity_k(3, Nband_Fermi_Level))
allocate( velocity_bar(3, Nk_start:Nk_end, Nband_Fermi_Level))
allocate( velocity_bar_k(3))
Ek= 0d0
Enk= 0d0
velocity= 0d0
velocity_k= 0d0
velocity_bar= 0d0
velocity_bar_k= 0d0
time_start= 0d0
time_end= 0d0
do ik= Nk_start, Nk_end
if (cpuid.eq.0.and. mod(ik, 100).eq.0) &
write(stdout, '(a, i18, " /", i18, a, f10.3, "s")') 'ik/NK', &
ik-Nk_start,Nk_current, 'time left', &
(Nk_current+Nk_start-ik)*(time_end-time_start)/num_cpu
call now(time_start)
k= KCube3D_total%k_direct(:, ik)
call velocity_calc(Nband_Fermi_Level, bands_fermi_level, k, velocity_k, Ek)
velocity(:, ik, :)= velocity_k
Enk(ik, :)= Ek
call now(time_end)
enddo
!> we get the kpath by Btau_final=-exponent_max*BTauMax, but we only use half of them
!> means that we can reach the accuracy as to exp(-exponent_max)
exponent_max= 30
!> exclude all kpoints with zero velocity x B and large energy away from Fermi level
do iband= 1, Nband_Fermi_Level
!> first check howmany k points left
it= 0
do ik= Nk_start, Nk_end
!> check whether v x B=0, which means the magnetic field is parallel with velocity
v_t= velocity(:, ik, iband)
!vcrossB(1)= -v_t(2)*Bdirection(3)+ v_t(3)*Bdirection(2)
!vcrossB(2)= -v_t(3)*Bdirection(1)+ v_t(1)*Bdirection(3)
!vcrossB(3)= -v_t(1)*Bdirection(2)+ v_t(2)*Bdirection(1)
!if (abs(Enk(ik, iband))<0.05d0.and.dsqrt(sum((abs(vcrossB)**2)))>eps3) then
if (abs(Enk(ik, iband))/eV2Hartree<EF_broadening) then
it = it+ 1
endif
enddo ! ik
KCube3D_left(iband)%Nk_current= it
if (it>0) then
allocate(KCube3D_left(iband)%ik_array(it))
allocate(KCube3D_left(iband)%Ek_local(it))
allocate(KCube3D_left(iband)%vx_local(it))
allocate(KCube3D_left(iband)%vy_local(it))
allocate(KCube3D_left(iband)%vz_local(it))
allocate(KCube3D_left(iband)%weight_k_local(it))
else
allocate(KCube3D_left(iband)%ik_array(1)) !> only useful for mpi_allgatherv
allocate(KCube3D_left(iband)%Ek_local(1))
allocate(KCube3D_left(iband)%vx_local(1))
allocate(KCube3D_left(iband)%vy_local(1))
allocate(KCube3D_left(iband)%vz_local(1))
allocate(KCube3D_left(iband)%weight_k_local(1))
endif
KCube3D_left(iband)%ik_array= 0
KCube3D_left(iband)%Ek_local= 0d0
KCube3D_left(iband)%vx_local= 0d0
KCube3D_left(iband)%vy_local= 0d0
KCube3D_left(iband)%vz_local= 0d0
KCube3D_left(iband)%weight_k_local= 0d0
it= 0
do ik= Nk_start, Nk_end
!> check whether v x B=0, which means the magnetic field is parallel with velocity
v_t= velocity(:, ik, iband)
!vcrossB(1)= -v_t(2)*Bdirection(3)+ v_t(3)*Bdirection(2)
!vcrossB(2)= -v_t(3)*Bdirection(1)+ v_t(1)*Bdirection(3)
!vcrossB(3)= -v_t(1)*Bdirection(2)+ v_t(2)*Bdirection(1)
!if (abs(Enk(ik, iband))<0.05d0.and.dsqrt(sum((abs(vcrossB)**2)))>eps3) then
if (abs(Enk(ik, iband))/eV2Hartree<EF_broadening) then
it = it+ 1
KCube3D_left(iband)%weight_k_local(it) = KCube3D_symm%weight_k(ik)
KCube3D_left(iband)%ik_array(it) = KCube3D_symm%ik_array_symm(ik)
KCube3D_left(iband)%Ek_local(it) = Enk(ik, iband)
KCube3D_left(iband)%vx_local(it) = v_t(1)
KCube3D_left(iband)%vy_local(it) = v_t(2)
KCube3D_left(iband)%vz_local(it) = v_t(3)
endif
enddo ! ik
enddo ! iband
!> try to get the total number of k points left for each band
do iband=1, Nband_Fermi_Level
#if defined (MPI)
call mpi_allreduce(KCube3D_left(iband)%Nk_current,KCube3D_left(iband)%Nk_total,1,&
mpi_in,mpi_sum,mpi_cmw,ierr)
#else
KCube3D_left(iband)%Nk_total= KCube3D_left(iband)%Nk_current
#endif
enddo
!> gather the number of k points left into a array info
allocate(info(num_cpu, Nband_Fermi_Level))
info= -1
do iband= 1, Nband_Fermi_Level
nrecevs= KCube3D_left(iband)%Nk_current
if (nrecevs<0) nrecevs= 0
#if defined (MPI)
call mpi_allgather(nrecevs, 1, mpi_in, info(:, iband), 1, mpi_in, mpi_cmw, ierr)
#else
info(1, iband)= KCube3D_left(iband)%Nk_current
#endif
enddo
!> An array for mpi_allgatherv
allocate(Displs(num_cpu+1, Nband_Fermi_Level))
Displs= 0
do iband=1, Nband_Fermi_Level
Displs(1, iband)=0
do i=2, num_cpu+1
Displs(i, iband)=Displs(i-1, iband)+ info(i-1, iband)
enddo
enddo ! iband
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
!> put all the kpoints left together
do iband=1, Nband_Fermi_Level
allocate(KCube3D_left(iband)%IKleft_array(KCube3D_left(iband)%Nk_total))
allocate(KCube3D_left(iband)%Ek_total(KCube3D_left(iband)%Nk_total))
allocate(KCube3D_left(iband)%vx_total(KCube3D_left(iband)%Nk_total))
allocate(KCube3D_left(iband)%vy_total(KCube3D_left(iband)%Nk_total))
allocate(KCube3D_left(iband)%vz_total(KCube3D_left(iband)%Nk_total))
allocate(KCube3D_left(iband)%weight_k(KCube3D_left(iband)%Nk_total))
KCube3D_left(iband)%IKleft_array = 0
KCube3D_left(iband)%Ek_total= 0d0
KCube3D_left(iband)%vx_total= 0d0
KCube3D_left(iband)%vy_total= 0d0
KCube3D_left(iband)%vz_total= 0d0
KCube3D_left(iband)%weight_k= 0d0
enddo ! iband
!> gather Enk and velocity
#if defined (MPI)
do iband=1, Nband_Fermi_Level
!nrecevs= KCube3D_left(iband)%Nk_end-KCube3D_left(iband)%Nk_start+ 1
nrecevs= KCube3D_left(iband)%Nk_current
if (nrecevs<0) nrecevs= 0
call mpi_allgatherv(KCube3D_left(iband)%ik_array, nrecevs, &
mpi_in, KCube3D_left(iband)%IKleft_array, &
info(:, iband), Displs(:, iband), mpi_in, mpi_cmw, ierr)
call mpi_allgatherv(KCube3D_left(iband)%Ek_local, nrecevs, &
mpi_dp, KCube3D_left(iband)%Ek_total, &
info(:, iband), Displs(:, iband), mpi_dp, mpi_cmw, ierr)
call mpi_allgatherv(KCube3D_left(iband)%vx_local, nrecevs, &
mpi_dp, KCube3D_left(iband)%vx_total, &
info(:, iband), Displs(:, iband), mpi_dp, mpi_cmw, ierr)
call mpi_allgatherv(KCube3D_left(iband)%vy_local, nrecevs, &
mpi_dp, KCube3D_left(iband)%vy_total, &
info(:, iband), Displs(:, iband), mpi_dp, mpi_cmw, ierr)
call mpi_allgatherv(KCube3D_left(iband)%vz_local, nrecevs, &
mpi_dp, KCube3D_left(iband)%vz_total, &
info(:, iband), Displs(:, iband), mpi_dp, mpi_cmw, ierr)
call mpi_allgatherv(KCube3D_left(iband)%weight_k_local, nrecevs, &
mpi_dp, KCube3D_left(iband)%weight_k, &
info(:, iband), Displs(:, iband), mpi_dp, mpi_cmw, ierr)
enddo ! iband
#else
do iband=1, Nband_Fermi_Level
KCube3D_left(iband)%IKleft_array= KCube3D_left(iband)%ik_array
KCube3D_left(iband)%Ek_total= KCube3D_left(iband)%Ek_local
KCube3D_left(iband)%vx_total= KCube3D_left(iband)%vx_local
KCube3D_left(iband)%vy_total= KCube3D_left(iband)%vy_local
KCube3D_left(iband)%vz_total= KCube3D_left(iband)%vz_local
KCube3D_left(iband)%weight_k= KCube3D_left(iband)%weight_k_local
enddo ! iband
#endif
!> redistribute all those k points into different cpus
if (cpuid.eq.0) then
write(stdout, '(a)')' '
write(stdout, '(a, i10, a)')' There are ', KCube3D_total%Nk_total, ' k points generated by the input files'
do iband= 1, Nband_Fermi_Level
write(stdout, '(a, i10, a, i10)')' However there are only ', KCube3D_left(iband)%Nk_total, &
' k points contribute to the conductance calculations for band ', bands_fermi_level(iband)
enddo
endif
!> redistribute the left kpoints into different processors
!> for different bands, the number of left kpoints is different.
do iband= 1, Nband_Fermi_Level
knv3_left= KCube3D_left(iband)%Nk_total
knv3_left_mod= mod(knv3_left, num_cpu)
if (knv3_left_mod==0) then !> perfect divided
KCube3D_left(iband)%Nk_current= knv3_left/num_cpu
KCube3D_left(iband)%Nk_start=1+ knv3_left*cpuid/num_cpu
KCube3D_left(iband)%Nk_end =(1+cpuid)*knv3_left/num_cpu
else if (knv3_left/num_cpu==0) then !> Number of MPI threads is large than knv3_left
KCube3D_left(iband)%Nk_current= 1 !> one k piont per MPI thread
KCube3D_left(iband)%Nk_start= cpuid+ 1 !> one k piont per MPI thread
KCube3D_left(iband)%Nk_end = cpuid+ 1
if (cpuid+1 > knv3_left) then
KCube3D_left(iband)%Nk_start= 1
KCube3D_left(iband)%Nk_end = 0
endif
else
KCube3D_left(iband)%Nk_current= knv3_left/num_cpu+ 1
if (cpuid< knv3_left_mod) then
KCube3D_left(iband)%Nk_start= 1+ cpuid*KCube3D_left(iband)%Nk_current
KCube3D_left(iband)%Nk_end = (1+cpuid)*KCube3D_left(iband)%Nk_current
else
KCube3D_left(iband)%Nk_start= knv3_left_mod*KCube3D_left(iband)%Nk_current+ &
(cpuid-knv3_left_mod)*(KCube3D_left(iband)%Nk_current-1)+1
KCube3D_left(iband)%Nk_end = knv3_left_mod*KCube3D_left(iband)%Nk_current+ &
(cpuid-knv3_left_mod+1)*(KCube3D_left(iband)%Nk_current-1)
endif
endif
allocate(KCube3D_left(iband)%k_direct(3, KCube3D_left(iband)%Nk_start:KCube3D_left(iband)%Nk_end))
do ik= KCube3D_left(iband)%Nk_start, KCube3D_left(iband)%Nk_end
i= KCube3D_left(iband)%IKleft_array(ik)
ik1= (i-1)/(Nk2*Nk3)+1
ik2= ((i-1-(ik1-1)*Nk2*Nk3)/Nk3)+1
ik3= (i-(ik2-1)*Nk3- (ik1-1)*Nk2*Nk3)
KCube3D_left(iband)%k_direct(1, ik)= (ik1-1)/dble(Nk1)
KCube3D_left(iband)%k_direct(2, ik)= (ik2-1)/dble(Nk2)
KCube3D_left(iband)%k_direct(3, ik)= (ik3-1)/dble(Nk3)
enddo ! ik
!> allocate array to store the conductivity for eack k point
allocate( sigma_iband_k(iband)%sigma_ohe_tensor_k(9, NBTau, OmegaNum, NumT), stat= ierr)
if (ierr>0) stop ' Error : no enough memory'
allocate( sigma_iband_k(iband)%sigma_ohe_tensor_k_mpi(9, NBTau, OmegaNum, NumT), stat= ierr)
if (ierr>0) stop ' Error : no enough memory'
! allocate( sigma_iband_k(iband)%sigma_ohe_tensor_kz(9, NBTau, OmegaNum, NumT, Nk3))
allocate( sigma_iband_k(iband)%time_cost(KCube3D_left(iband)%Nk_total))
allocate( sigma_iband_k(iband)%time_cost_mpi(KCube3D_left(iband)%Nk_total))
sigma_iband_k(iband)%time_cost= 0d0
sigma_iband_k(iband)%time_cost_mpi= 0d0
sigma_iband_k(iband)%sigma_ohe_tensor_k= 0d0
! sigma_iband_k(iband)%sigma_ohe_tensor_kz= 0d0
sigma_iband_k(iband)%sigma_ohe_tensor_k_mpi= 0d0
enddo ! iband=1, Nband_Fermi_Level
!> gather the number of receive buffs left into a array info
info= -1
do iband= 1, Nband_Fermi_Level
nrecevs= (KCube3D_left(iband)%Nk_end-KCube3D_left(iband)%Nk_start+1)*9*NBTau*OmegaNum*NumT
if (nrecevs<0) nrecevs=0
#if defined (MPI)
call mpi_allgather(nrecevs, 1, mpi_in, info(:, iband), 1, mpi_in, mpi_cmw, ierr)
#else
info(1, iband)= nrecevs
#endif
enddo
if (cpuid.eq.0) write(stdout, '(a, 10i8)') ' '
do iband=1, Nband_Fermi_Level
if (cpuid.eq.0) write(stdout, '(a, i7)')'>> Number of k points at different CPUs at band', iband
if (cpuid.eq.0) write(stdout, '(10i8)')info(:, iband)/9/NBTau/OmegaNum/NumT
enddo
if (cpuid.eq.0) write(stdout, '(a, 10i8)') ' '
!> An array for mpi_allgatherv
Displs= 0
do iband=1, Nband_Fermi_Level
Displs(1, iband)=0
do i=2, num_cpu+1
Displs(i, iband)=Displs(i-1, iband)+ info(i-1, iband)
enddo
enddo ! iband
do iband=1, Nband_Fermi_Level
allocate(klist_iband(iband)%klist_rkfs(3, NSlice_Btau))
allocate(klist_iband(iband)%velocity_k(3, NSlice_Btau))
klist_iband(iband)%klist_rkfs= 0d0
klist_iband(iband)%velocity_k= 0d0
enddo
!> a temp array used in RKFS
allocate(kout(3, NSlice_Btau))
kout= 0d0
!> allocate array for sigmafileindex and rhofileindex
allocate(sigmafileindex(OmegaNum))
allocate(rhofileindex(OmegaNum))
sigmafileindex= (/(outfileindex + ie, ie = 1, OmegaNum)/)
rhofileindex= (/(outfileindex + ie +OmegaNum, ie = 1, OmegaNum)/)
outfileindex = outfileindex + 2*OmegaNum
if (cpuid.eq.0) then
do ie = 1, OmegaNum
write(muname, '(f12.2)')mu_array(ie)/eV2Hartree
!> open file for conductivity at mu and write header
write(sigmafilename, '(3a)')'sigma_bands_mu_',trim(adjustl(muname)),'eV.dat'
open(unit=sigmafileindex(ie), file=sigmafilename)
write(sigmafileindex(ie), '(a)')'# Conductivity tensor/tau (in unit of (\Omega*m*s)^-1) for every contributing band'
write(sigmafileindex(ie), '(a,100I5)')'# SELECTEDBANDS = ', bands_fermi_level(:)
write(sigmafileindex(ie), '(a,1000f8.3)')'# Tlist = ', KBT_array(:)/8.6173324E-5/eV2Hartree
write(sigmafileindex(ie), '(a, 1000E16.6)') '# Btaulist = ', BTau_array(:)*Magneticfluxdensity_atomic/Relaxation_Time_Tau
!> open file for resistivity at mu and write header
write(sigmafilename, '(3a)')'rho_bands_mu_', trim(adjustl(muname)),'eV.dat'
open(unit=rhofileindex(ie), file=sigmafilename)
write(rhofileindex(ie), '(a)')'# Resistivity \tau*\rho (in unit of \Omega*m*s) for every contributing band '
write(rhofileindex(ie), '(a,100I5)')'# SELECTEDBANDS = ', bands_fermi_level(:)
write(rhofileindex(ie), '(a,1000f8.3)')'# Tlist = ', KBT_array(:)/8.6173324E-5/eV2Hartree
write(rhofileindex(ie), '(a, 1000E16.6)') '# Btaulist = ', BTau_array(:)*Magneticfluxdensity_atomic/Relaxation_Time_Tau
enddo ! ie
endif ! cpuid.eq.0
!> now we turn to use Runge-Kutta method to get all the kpoints from (0, BTauMax)
!> and we calculate the conductivity/Tau over different bands and different k points
time_start= 0d0
time_end = 0d0
do iband= 1, Nband_Fermi_Level
call now(time_start0)
!> dim=(Nk_start: Nk_end, NumT, OmegaNum))
allocate(klist_iband(iband)%minusdfde(OmegaNum, NumT))
call cal_sigma_iband_k
#if defined (MPI)
call mpi_allreduce(sigma_iband_k(iband)%time_cost_mpi, sigma_iband_k(iband)%time_cost, &
size(sigma_iband_k(iband)%time_cost), &
mpi_dp,mpi_sum,mpi_cmw,ierr)
call mpi_allreduce(sigma_iband_k(iband)%sigma_ohe_tensor_k_mpi, sigma_iband_k(iband)%sigma_ohe_tensor_k, &
size(sigma_iband_k(iband)%sigma_ohe_tensor_k), mpi_dp,mpi_sum,mpi_cmw,ierr)
if (ierr>0) then
write(stdout, *)'>>>Error happends in mpi_allgatherv at cpuid', cpuid, ' ierr=', ierr
stop
endif
#else
sigma_iband_k(iband)%time_cost= sigma_iband_k(iband)%time_cost_mpi
sigma_iband_k(iband)%sigma_ohe_tensor_k= sigma_iband_k(iband)%sigma_ohe_tensor_k_mpi
#endif
!call mpi_barrier(mpi_cmw, ierr)
!call now(time_2)
!if (cpuid.eq.0) write(stdout, *) 'after mpi_allgatherv in sigma_ohe_calc', time_2-time_1, 's'
if (cpuid.eq.0) then
write(stdout, '(a)')' '
write(stdout, '(a, i10)')'>> Time cost for each k point at iband ', iband
write(stdout, '(10f16.2)')(sigma_iband_k(iband)%time_cost(ik), ik= 1, KCube3D_left(iband)%Nk_total)
write(stdout, '(a)')' '
endif
!> sum over all the kpoints to get the band dependent conductivity/tau
do ikt=1, NumT
do ie=1, OmegaNum
do ibtau=1, NBTau
do i=1, 9
sigma_ohe_tensor(i, ibtau, ie, ikt, iband)= &
sigma_iband_k(iband)%sigma_ohe_tensor_k(i, ibtau, ie, ikt)
!print *, i, ibtau, ie, ikt, sigma_ohe_tensor(i, ibtau, ie, ikt, iband)
!pause
enddo ! i
enddo ! ibtau
enddo ! ie
enddo ! ikt
call now(time_end0)
if (cpuid.eq.0) write(stdout, '(a, i6, a, f16.2, a)')'>> Time cost for calculate MR at iband=', iband, &
'is ', time_end0- time_start0, ' s'
!> In the end, we start to care about the units of the conductivity/tau
!> change from the Hatree Atomic units to SI units
!> the conductivity/tau is in units of Ohm^-1*m^-1*s^-1
!>
coeff= Echarge**2/hbar/Bohr_radius/Origin_cell%CellVolume/kCubeVolume*Origin_cell%ReciprocalCellVolume
coeff= coeff/Time_atomic
sigma_ohe_tensor(:, :, :, :, iband)= sigma_ohe_tensor(:, :, :, :, iband)*coeff
if (cpuid.eq.0) then
do ie = 1, OmegaNum
!> open file for conductivity at mu and write data
write(sigmafileindex(ie), '(2a, i5)')'# ',' iband = ', bands_fermi_level(iband)
write(sigmafileindex(ie),'(a)') ''
do ikt = 1, NumT
KBT= KBT_array(ikt)/8.6173324E-5/eV2Hartree
write(sigmafileindex(ie), '(2a, f16.4, a)') '#', ' T = ', KBT, ' K'
write(sigmafileindex(ie), '("# Column", i5, 100i16)')(i, i=1, 10)
write(sigmafileindex(ie), '("#",19a16)')'BTau (T.ps)', 'xx', 'xy', 'xz', 'yx', 'yy', 'yz', 'zx', 'zy','zz'
!> write out the conductivity/tau into file
do i=1, NBTau
write(sigmafileindex(ie), '(19E16.6)') &
BTau_array(i)*Magneticfluxdensity_atomic/Relaxation_Time_Tau, &
sigma_ohe_tensor(:, i, ie, ikt, iband)
enddo ! i, NBTau
write(sigmafileindex(ie),'(a)') ''
enddo
write(sigmafileindex(ie),'(a)') ''
!> open file for resistivity at mu and write data
write(rhofileindex(ie), '(2a, i5)')'# ',' iband = ', bands_fermi_level(iband)
write(rhofileindex(ie),'(a)') ''
do ikt = 1, NumT
KBT= KBT_array(ikt)/8.6173324E-5/eV2Hartree
write(rhofileindex(ie), '(2a, f16.4, a)') '#', ' T = ', KBT, ' K'
write(rhofileindex(ie), '("# Column", i5, 100i16)')(i, i=1, 10)
write(rhofileindex(ie), '("#",19a16)')'BTau (T.ps)', 'xx', 'xy', 'xz', 'yx', 'yy', 'yz', 'zx', 'zy','zz'
!> write out the inverse of conductivity/tau into file
do i=1, NBTau
rho(1, 1:3)=sigma_ohe_tensor(1:3, i, ie, ikt, iband)
rho(2, 1:3)=sigma_ohe_tensor(4:6, i, ie, ikt, iband)
rho(3, 1:3)=sigma_ohe_tensor(7:9, i, ie, ikt, iband)
if (abs(det3(rho))>eps6) then
call inv_r(3, rho)
write(rhofileindex(ie), '(19E16.6)')&
BTau_array(i)*Magneticfluxdensity_atomic/Relaxation_Time_Tau, &
rho(1, 1:3), rho(2, 1:3), rho(3, 1:3)
else
write(rhofileindex(ie), '(a)')'# error: sigma is zero since no k points contribute to the calculations of MR'
endif
enddo ! i, NBTau
write(rhofileindex(ie),'(a)')''
enddo
write(rhofileindex(ie),'(a)')''
enddo ! ie, OmegaNum
endif
!---------------------------------------------------------------------------------------------------------------
!> This part is to calculate and store sigma_kz and was commentted out by pihanqi on Apr. 29. 2023
! !calculate sigma_kz
! do ik3=1, Nk3
! do ik= 1, KCube3D_left(iband)%Nk_total
! i= KCube3D_left(iband)%IKleft_array(ik)
! ik1= (i-1)/(Nk2*Nk3)+1
! ik2= ((i-1-(ik1-1)*Nk2*Nk3)/Nk3)+1
! if (ik3== (i-(ik2-1)*Nk3- (ik1-1)*Nk2*Nk3)) then
! do ikt=1, NumT
! do ie=1, OmegaNum
! do ibtau=1, NBTau
! do ix=1, 9
! sigma_iband_k(iband)%sigma_ohe_tensor_kz(ix, ibtau, ie, ikt, ik3)= &
! sigma_iband_k(iband)%sigma_ohe_tensor_kz(ix, ibtau, ie, ikt, ik3)+ &
! sigma_iband_k(iband)%sigma_ohe_tensor_k (ix, ibtau, ie, ikt, ik )
! enddo ! ix
! enddo ! ibtau
! enddo ! ie
! enddo ! ikt
! endif
! enddo ! ik
! enddo ! ik3
! sigma_iband_k(iband)%sigma_ohe_tensor_kz= &
! sigma_iband_k(iband)%sigma_ohe_tensor_kz*coeff
!allocate(myfileindex(Nband_Fermi_Level))
!!> file index for different bands
!do iband=1, Nband_Fermi_Level
! outfileindex= outfileindex+ 1
! myfileindex(iband)= outfileindex
!enddo
! do ikt=1, NumT
! do ie=1, OmegaNum
! outfileindex= outfileindex+ 1
! KBT= KBT_array(ikt)/8.6173324E-5/eV2Hartree
! write(tname, '(f12.2)')KBT
! write(muname, '(f12.2)')mu_array(ie)/eV2Hartree
! if (cpuid.eq.0) then
! write(bandname, '(i10)')bands_fermi_level(iband)
! write(sigmafilename, '(7a)')'sigma_kz_band_', trim(adjustl(bandname)),'_mu_',&
! trim(adjustl(muname)),'eV_T_', trim(adjustl(tname)), 'K.dat'
! open(unit=outfileindex, file=sigmafilename)
! write(outfileindex, '(a20, i5, 2(a, f16.4, a))')'# sigma_k at band ', &
! bands_fermi_level(iband), ' temperature at ', KBT, ' K', ' chemical potential at ', mu_array(ie), ' eV'
! write(outfileindex, '("#",a12,20a16)')'ik3', 'xx', 'xy', 'xz', 'yx', 'yy', 'yz', 'zx', 'zy','zz'
! endif
! if (cpuid.eq.0) then
! do ibtau=1, NBTau
! !> write out the conductivity/tau into file
! write(outfileindex, '(a6, f16.8, a12, f16.8)')'# Btau', &
! BTau_array(ibtau)*Magneticfluxdensity_atomic/Relaxation_Time_Tau, &
! ' omega*tau', BTau_array(ibtau)*Magneticfluxdensity_atomic/Relaxation_Time_Tau*0.175874356d0
! write(outfileindex, '("#",a12,20a16)')'ik3', 'xx', 'xy', 'xz', 'yx', 'yy', 'yz', 'zx', 'zy','zz'
! do ik3=1, Nk3
! write(outfileindex, '(i6,20E16.6)')ik3, sigma_iband_k(iband)%sigma_ohe_tensor_kz(:, ibtau, ie, ikt, ik3)
! enddo
! enddo ! ibtau
! endif ! cpuid=0
! if (cpuid.eq.0) then
! close(outfileindex)
! endif ! cpuid=0
!---------------------------------------------------------------------------------------------------------------
enddo ! iband =1, Nband_Fermi_Level
if (cpuid.eq.0) then
do ie = 1, OmegaNum
!> close file for conductivity at mu
close(sigmafileindex(ie))
!> close file for resistivity at mu
close(rhofileindex(ie))
enddo ! ie
endif ! cpuid.eq.0
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
!> write data into file, modified by pihanqi on Apr. 29. 2023
!> write out the conductivity/tau and resistivity*tau with the same relaxation time for every contributing bands
if (cpuid.eq.0) then
do ie=1, OmegaNum
write(muname, '(f12.2)')mu_array(ie)/eV2Hartree
!> write out the total conductivity with the same relaxation time for all bands
outfileindex= outfileindex+ 1
write(sigmafilename, '(3a)')'sigma_total_mu_',trim(adjustl(muname)),'eV.dat'
open(unit=outfileindex, file=sigmafilename)
write(outfileindex, '(a)')'# \sigma/\tau with unit (Ohm*m*s)^-1 is the summation of all bands \sum_n\sigma_n/\tau_n '
write(outfileindex, '(a)')'# relaxation time \tau_n=\tau is the same for all bands. '
write(outfileindex, '(a,2I6)')'# NumT NumBtau = ', NumT, NBTau
write(outfileindex, '(a,1000f8.3)')'# Tlist = ', KBT_array(:)/8.6173324E-5/eV2Hartree
do ikt = 1, NumT
KBT= KBT_array(ikt)/8.6173324E-5/eV2Hartree
write(outfileindex, '(2a, f16.4, a)')'# ',' T = ', KBT, ' K '
write(outfileindex, '("# Column", i5, 100i16)')(i, i=1, 10)
write(outfileindex, '("#",19a16)')'BTau (T.ps)', 'xx', 'xy', 'xz', 'yx', 'yy', 'yz', 'zx', 'zy','zz'
!> write out the inverse of conductivity/tau into file
!> the name of rho is meaningless here, just a temp variable
do i=1, NBTau
rho(1, 1:3)=sum(sigma_ohe_tensor(1:3, i, ie, ikt, :), dim=2)
rho(2, 1:3)=sum(sigma_ohe_tensor(4:6, i, ie, ikt, :), dim=2)
rho(3, 1:3)=sum(sigma_ohe_tensor(7:9, i, ie, ikt, :), dim=2)
write(outfileindex, '(19E16.6)')&
BTau_array(i)*Magneticfluxdensity_atomic/Relaxation_Time_Tau, &
rho(1, 1:3), rho(2, 1:3), rho(3, 1:3)
enddo ! i, NBTau
write(outfileindex, '(a, 10i8)')''
enddo
close(outfileindex)
!> write out the total resistivity with the same relaxation time for all bands
outfileindex= outfileindex+ 1
write(sigmafilename, '(3a)')'rho_total_mu_',trim(adjustl(muname)),'eV.dat'
open(unit=outfileindex, file=sigmafilename)
write(outfileindex, '(a)')'# \tau*\rho with unit (Ohm*m*s) is the inverse of Conductivity tensor \sum_n\sigma_n/\tau '
write(outfileindex, '(a)')'# relaxation time \tau_n=\tau is the same for all bands. '
write(outfileindex, '(a,2I6)')'# NumT NumBtau = ', NumT, NBTau
write(outfileindex, '(a,1000f8.3)')'# Tlist = ', KBT_array(:)/8.6173324E-5/eV2Hartree
do ikt =1, NumT
KBT= KBT_array(ikt)/8.6173324E-5/eV2Hartree
write(outfileindex, '(2a, f16.4, a)')'# ', ' T = ', KBT, ' K '
write(outfileindex, '("# Column", i5, 100i16)')(i, i=1, 10)
write(outfileindex, '("#",19a16)')'BTau (T.ps)', 'xx', 'xy', 'xz', 'yx', 'yy', 'yz', 'zx', 'zy','zz'
!> write out the inverse of conductivity/tau into file
do i=1, NBTau
rho(1, 1:3)=sum(sigma_ohe_tensor(1:3, i, ie, ikt, :), dim=2)
rho(2, 1:3)=sum(sigma_ohe_tensor(4:6, i, ie, ikt, :), dim=2)
rho(3, 1:3)=sum(sigma_ohe_tensor(7:9, i, ie, ikt, :), dim=2)
if (abs(det3(rho))>eps6) then
call inv_r(3, rho)
write(outfileindex, '(19E16.6)')&
BTau_array(i)*Magneticfluxdensity_atomic/Relaxation_Time_Tau, &
rho(1, 1:3), rho(2, 1:3), rho(3, 1:3)
else
write(outfileindex, '(a)')'# error: sigma is zero since no k points contribute to the calculations of MR'
endif
enddo ! i, NBTau
write(outfileindex, '(a, 10i8)')''
enddo
close(outfileindex)
enddo ! ie=1, OmegaNum
endif ! cpuid=0
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
contains
subroutine cal_sigma_iband_k
do ik= KCube3D_left(iband)%Nk_start, KCube3D_left(iband)%Nk_end
if (cpuid.eq.0) &
write(stdout, '(a, i8, a, i18, " /", i18, a, f10.3, "s", a, f10.3, "s")') &
'In sigma_OHE iband', iband, ' ik/NK', &
ik-KCube3D_left(iband)%Nk_start+1,KCube3D_left(iband)%Nk_current, &
' time cost', time_end-time_start, &
' time left', &
(KCube3D_left(iband)%Nk_current+KCube3D_left(iband)%Nk_start-ik)*(time_end-time_start)
call now(time_start)
EE= KCube3D_left(iband)%Ek_total(ik)
v_k(1)= KCube3D_left(iband)%vx_total(ik)
v_k(2)= KCube3D_left(iband)%vy_total(ik)
v_k(3)= KCube3D_left(iband)%vz_total(ik)
!> calculate df/de for each k point and each band
do ikt=1, NumT
KBT= KBT_array(ikt)
do ie=1, OmegaNum
mu= mu_array(ie)
call minusdfde_calc_single(EE, KBT, mu, minusdfde)
klist_iband(iband)%minusdfde(ie, ikt)= minusdfde
enddo ! ie
enddo ! ikt
!> start to get the evolution of k points under magnetic field using Runge-Kutta method
k_start= KCube3D_left(iband)%k_direct(:, ik)
kout= 0d0
Btau_start= 0d0
!> we get the kpath by Btau_final=-30*BTauMax, but we only use half of them
Btau_final= -exponent_max*BTauMax !< -15 means that we can reach the accuracy as to exp(-15d0)
!Btau_final= -10d0*BTauMax !< test for arXiv:2201.03292 (2022)
!> Runge-Kutta only applied with BTauMax>0
!> if the magnetic field is zero.
if (BTauMax>eps3) then
NSlice_Btau_inuse= NSlice_Btau
call RKF45_pack(magnetic_field, bands_fermi_level(iband), &
NSlice_Btau_inuse, k_start, Btau_start, Btau_final, kout, icycle, fail)
else
icycle= 1
do ibtau=1, NSlice_Btau
kout(:, ibtau)= k_start(:)
enddo
NSlice_Btau_inuse = NSlice_Btau
endif
if (NSlice_Btau_inuse==1) then
write(stdout, '(a, i6, a, i4, a, i6, a, 3f12.6)')&
'>>> NSlice_Btau_inuse=1 at cpuid=', cpuid, ' iband=', iband, ' ik', ik, ' k', k_start
endif
!> we omit the
if (fail) then
write(stdout, '(a, i6, a, i4, a, i6, a, 3f12.6)')&
'>>> Runge-Kutta integration fails at cpuid=', cpuid, ' iband=', iband, ' ik', ik, ' k', k_start
write(stdout, *)' '
cycle
endif
if (NSlice_Btau_inuse==1)then ! vxB=0
call velocity_calc_iband(bands_fermi_level(iband), k_start, v_t)
do ibtau=1, NSlice_Btau
kout(:, ibtau)= k_start(:)
klist_iband(iband)%velocity_k(:, ibtau)= v_t
enddo
NSlice_Btau_inuse = NSlice_Btau
else
do it= 1, icycle
k= kout(:, it)
call velocity_calc_iband(bands_fermi_level(iband), k, v_t)
klist_iband(iband)%velocity_k(:, it)= v_t
enddo ! integrate over time step
!> periodic kpath in the BZ can be reused
do i=2, NSlice_Btau_inuse/icycle
do j=1, icycle
klist_iband(iband)%velocity_k(:, j+(i-1)*icycle)= klist_iband(iband)%velocity_k(:, j)
enddo
enddo
do i=(NSlice_Btau_inuse/icycle)*icycle+1, NSlice_Btau
klist_iband(iband)%velocity_k(:, i)= klist_iband(iband)%velocity_k(:, i-(NSlice_Btau_inuse/icycle)*icycle)
enddo
endif
!> calculate the conductivity/tau
do ikt=1, NumT
KBT= KBT_array(ikt)
do ie=1, OmegaNum
mu= mu_array(ie)
minusdfde= klist_iband(iband)%minusdfde(ie, ikt)
do ibtau=1, NBTau
BTau= BTau_array(ibtau)
if (NBTau==1)then
NSlice_Btau_local= 2
else
NSlice_Btau_local= (ibtau-1)*NSlice_Btau_inuse/(NBTau-1)/2
if (NSlice_Btau_local==0)then
NSlice_Btau_local= 2
else
DeltaBtau= exponent_max/2d0/NSlice_Btau_local
endif
endif
!do ishift=0, NSlice_Btau_local-1
do ishift=0, 0
!> here, velocity is in the cartesian coordinate
!> The core of Chamber formular is to get the average of velocity
!> in the relaxation time approximation
v_k= klist_iband(iband)%velocity_k(:, 1+ishift)
if (BTau>eps3.and.NSlice_Btau_local>2) then
velocity_bar_k= 0d0
do it=1, NSlice_Btau_local
velocity_bar_k= velocity_bar_k+ &
DeltaBtau*exp(-(it-1d0)*DeltaBtau)*klist_iband(iband)%velocity_k(:, it+ishift)
enddo
velocity_bar_k= velocity_bar_k &
- 0.5d0*DeltaBtau*(exp(-(NSlice_Btau_local-1d0)*DeltaBtau)&
* klist_iband(iband)%velocity_k(:, NSlice_Btau_local+ishift) &
+ klist_iband(iband)%velocity_k(:, 1+ishift))
else
velocity_bar_k= v_k
endif
!> calculate the conductivity now
ik_temp= ik - KCube3D_left(iband)%Nk_start+ 1
sigma_symm_t= 0d0
!> Apply point group operations to the velocities, and average them
do j1=1, 3
do j2=1, 3
do j=1, number_group_operators
!> sigma_xx
sigma_symm_t(1)= sigma_symm_t(1)+ &
pgop_cart_inverse(1, j1, j)* pgop_cart_inverse(1, j2, j) &
*v_k(j1)*velocity_bar_k(j2)
!> sigma_xy
sigma_symm_t(2)= sigma_symm_t(2)+ &
pgop_cart_inverse(1, j1, j)* pgop_cart_inverse(2, j2, j) &
*v_k(j1)*velocity_bar_k(j2)
!> sigma_xz
sigma_symm_t(3)= sigma_symm_t(3)+ &
pgop_cart_inverse(1, j1, j)* pgop_cart_inverse(3, j2, j) &
*v_k(j1)*velocity_bar_k(j2)
!> sigma_yx
sigma_symm_t(4)= sigma_symm_t(4)+ &
pgop_cart_inverse(2, j1, j)* pgop_cart_inverse(1, j2, j) &
*v_k(j1)*velocity_bar_k(j2)
!> sigma_yy
sigma_symm_t(5)= sigma_symm_t(5)+ &
pgop_cart_inverse(2, j1, j)* pgop_cart_inverse(2, j2, j) &
*v_k(j1)*velocity_bar_k(j2)
!> sigma_yz
sigma_symm_t(6)= sigma_symm_t(6)+ &
pgop_cart_inverse(2, j1, j)* pgop_cart_inverse(3, j2, j) &
*v_k(j1)*velocity_bar_k(j2)
!> sigma_zx
sigma_symm_t(7)= sigma_symm_t(7)+ &
pgop_cart_inverse(3, j1, j)* pgop_cart_inverse(1, j2, j) &