-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathfermisurface.f90
executable file
·1924 lines (1583 loc) · 66 KB
/
fermisurface.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 fermisurface3D
! This subroutine calculates 3D fermi surface in the 1st BZ using wannier TB method
!
! Copyright (c) 2010 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
integer :: ik, i, knv3, ikx, iky, ikz, ierr
integer :: nband_min, nband_max, nband_store
character(40) :: fsfile
real(dp) :: k(3)
real(dp) :: time_start, time_end
! Hamiltonian of bulk system
complex(Dp), allocatable :: Hamk_bulk(:, :)
real(dp) :: kxmin, kxmax, kymin, kymax, kzmin, kzmax
real(dp), allocatable :: W(:)
real(dp), allocatable :: eigval(:,:)
real(dp), allocatable :: eigval_mpi(:,:)
! only for output the FS3D.bxsf, we don't have to output all the bands,
! only consider the bands close to the Fermi level
if (SOC == 0) then
nband_min= Numoccupied- 5
nband_max= Numoccupied+ 6
else
nband_min= Numoccupied- 7
nband_max= Numoccupied+ 8
endif
if (nband_min< 1) then
nband_min= 1
endif
if (nband_max> Num_wann) then
nband_max= Num_wann
endif
if (nband_min>nband_max) then
nband_min= max(1, nband_max-4)
if (SOC>0) nband_min= max(1, nband_max-8)
endif
nband_store= nband_max- nband_min+ 1
if (cpuid.eq.0) then
write(stdout, *)">> In fermisurface3D"
write(stdout, '(a, i8, a, i8)')">> We are going to write out the bands from band ", nband_min, ' to', nband_max
endif
kxmin= 0.00d0/1d0
kxmax= 1.00d0/1d0
kymin= 0.00d0/1d0
kymax= 1.00d0/1d0
kzmin= 0.00d0/1d0
kzmax= 1.00d0/1d0
ik =0
knv3= nk1*nk2*nk3
allocate(W(Num_wann))
allocate(Hamk_bulk(Num_wann, Num_wann))
allocate(eigval(nband_store, knv3), stat=ierr)
if (ierr>0) stop 'not enough memory'
allocate(eigval_mpi(nband_store, knv3), stat=ierr)
if (ierr>0) stop 'not enough memory'
eigval_mpi= 0d0
eigval= 0d0
time_start= 0d0
time_end= 0d0
do ik= 1+cpuid, knv3, num_cpu
if (cpuid==0.and. mod(ik/num_cpu, 500)==0) &
write(stdout, *) '3DFS, ik ', ik, 'knv3',knv3, 'time left', &
(knv3-ik)*(time_end- time_start)/num_cpu, ' s'
call now(time_start)
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-1) &
+ K3D_vec2_cube*(iky-1)/dble(nk2-1) &
+ K3D_vec3_cube*(ikz-1)/dble(nk3-1)
! calculation bulk hamiltonian
Hamk_bulk= 0d0
!call ham_bulk_atomicgauge (k, Hamk_bulk)
call ham_bulk_latticegauge (k, Hamk_bulk)
call eigensystem_c( 'N', 'U', Num_wann, Hamk_bulk, W)
eigval_mpi(:, ik)= W(nband_min:nband_max)
call now(time_end)
enddo
#if defined (MPI)
call mpi_allreduce(eigval_mpi, eigval,size(eigval),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
if (ierr>0) then
print *, 'Something wrong in mpi_allreduce in fermisurface3D', ierr
stop
endif
#else
eigval= eigval_mpi
#endif
if (cpuid==0) then
write(stdout, *)'>> All processors finished their job for fermisurface3D'
endif
!> writeout eigenvalues data for each band into separate files for matlab isosurface function use.
do i=1, nband_store
outfileindex= outfileindex+ 1
if (cpuid==0) then
if (i>0 .and. i<10) then
write(fsfile, '(a, i1, a)')'FS3D_matlab_band_',i, '.txt'
else if (i>=10 .and. i<100) then
write(fsfile, '(a, i2, a)')'FS3D_matlab_band_',i, '.txt'
else if (i>=100 .and. i<1000) then
write(fsfile, '(a, i3, a)')'FS3D_matlab_band_',i, '.txt'
else if (i>=1000 .and. i<10000) then
write(fsfile, '(a, i4, a)')'FS3D_matlab_band_',i, '.txt'
endif
open(unit=outfileindex,FILE=fsfile,STATUS='UNKNOWN',FORM='FORMATTED')
write(outfileindex,'(a, i10)') '%BAND: ',i
write(outfileindex,'(a, 4i16)') '%Nk1, Nk2, Nk3, total', Nk1, Nk2, Nk3, knv3
do ik=1, knv3
write(outfileindex,'(f16.8)') eigval(i, ik)/eV2Hartree
enddo
close(outfileindex)
endif
enddo
outfileindex= outfileindex+ 1
if (cpuid==0)then
open(unit=outfileindex,FILE='FS3D.bxsf',STATUS='UNKNOWN',FORM='FORMATTED')
write(outfileindex,'(a)') ' BEGIN_INFO'
write(outfileindex,'(a)') ' #'
write(outfileindex,'(a)') ' # this is a Band-XCRYSDEN-Structure-File'
write(outfileindex,'(a)') ' # for Fermi Surface Visualisation'
write(outfileindex,'(a)') ' #'
write(outfileindex,'(a)') ' # Launch as: xcrysden --bxsf FS3D.bxsf'
write(outfileindex,'(a)') ' Fermi Energy: 0'
write(outfileindex,'(a)') ' END_INFO'
write(outfileindex,'(a)')
write(outfileindex,'(a)') ' BEGIN_BLOCK_BANDGRID_3D'
write(outfileindex,'(a)') 'from_wannier_code'
write(outfileindex,'(a)') ' BEGIN_BANDGRID_3D_fermi'
write(outfileindex,'(i10)') nband_store
write(outfileindex,'(3i10)') nk1, nk2, nk3
write(outfileindex,'(a)') '0.0 0.0 0.0'
write(outfileindex,'(3f16.8)') (Origin_cell%Kua(i)*Angstrom2atomic, i=1,3)
write(outfileindex,'(3f16.8)') (Origin_cell%Kub(i)*Angstrom2atomic, i=1,3)
write(outfileindex,'(3f16.8)') (Origin_cell%Kuc(i)*Angstrom2atomic, i=1,3)
do i=1,nband_store
write(outfileindex,'(a, i10)') 'BAND: ',i
do ik=1, knv3
write(outfileindex,'(E16.8)') eigval(i, ik)/eV2Hartree
enddo
enddo
write(outfileindex,'(a)') 'END_BANDGRID_3D'
write(outfileindex,'(a)') ' END_BLOCK_BANDGRID_3D'
close(outfileindex)
endif
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
deallocate(W)
deallocate(Hamk_bulk)
deallocate(eigval)
deallocate(eigval_mpi)
return
end subroutine fermisurface3D
subroutine orbitaltexture
! This subroutine calculates orbital texture in a k plane
! using wannier TB method
use wmpi
use para
implicit none
integer :: ia, ig
integer :: ik, i, j, i1, i2
integer :: knv3, nkx, nky
integer :: ierr
real(dp) :: kz
real(Dp) :: k(3)
real(Dp) :: dos_max
! Hamiltonian of bulk system
complex(Dp), allocatable :: Hamk_bulk(:, :)
real(dp) :: time_start, time_end
real(dp) :: kxmin, kxmax, kymin, kymax
real(dp) :: zmin, zmax
real(dp) :: kxmin_shape, kxmax_shape, kymin_shape, kymax_shape
real(dp), allocatable :: kxy(:,:)
real(dp), allocatable :: kxy_shape(:,:)
real(dp), allocatable :: kxy_plane(:,:)
!> atomic and k resolved spectral function
real(dp), allocatable :: dos(:, :)
real(dp), allocatable :: dos_mpi(:, :)
real(dp), allocatable :: dos_total(:)
real(dp), allocatable :: dos_selected(:,:)
real(dp), allocatable :: dos_unselected(:,:)
complex(dp), allocatable :: ones(:,:)
allocate(Hamk_bulk(Num_wann, Num_wann))
nkx= Nk
nky= Nk
allocate( kxy(3, nkx*nky))
allocate( kxy_shape(3, nkx*nky))
allocate( kxy_plane(3, nkx*nky))
kxy=0d0
kxy_shape=0d0
kxy_plane=0d0
ik =0
do i= 1, nkx
do j= 1, nky
ik =ik +1
kxy(:, ik)= K3D_start+ K3D_vec1*(i-1)/dble(nkx-1)+ K3D_vec2*(j-1)/dble(nky-1)
kxy_shape(:, ik)= kxy(1, ik)* Origin_cell%Kua+ kxy(2, ik)* Origin_cell%Kub+ kxy(3, ik)* Origin_cell%Kuc
call rotate_k3_to_kplane(kxy_shape(:, ik), kxy_plane(:, ik))
enddo
enddo
i1=1
i2=2
kxmin_shape=minval(kxy_shape(i1,:))
kxmax_shape=maxval(kxy_shape(i1,:))
kymin_shape=minval(kxy_shape(i2,:))
kymax_shape=maxval(kxy_shape(i2,:))
knv3= nkx*nky
allocate( dos (knv3, Num_wann))
allocate( dos_mpi(knv3, Num_wann))
allocate( dos_total(knv3))
allocate( dos_selected(knv3, NumberofSelectedOrbitals_groups))
allocate( dos_unselected(knv3, NumberofSelectedOrbitals_groups))
dos = 0d0
dos_mpi= 0d0
allocate(ones(Num_wann, Num_wann))
ones= 0d0
do i=1, Num_wann
ones(i, i)= 1d0
enddo
time_start= 0d0
time_end= 0d0
do ik= 1+cpuid, knv3, num_cpu
if (cpuid==0.and. mod(ik/num_cpu, 500)==0) &
write(stdout, *) 'FS_Plane, ik ', ik, 'knv3',knv3, 'time left', &
(knv3-ik)*(time_end- time_start)/num_cpu, ' s'
call now(time_start)
k = kxy(:, ik)
! calculation bulk hamiltonian
Hamk_bulk= 0d0
if (index(KPorTB, 'KP')/=0)then
call ham_bulk_kp(k, Hamk_bulk)
else
!> deal with phonon system
if (index(Particle,'phonon')/=0.and.LOTO_correction) then
call ham_bulk_LOTO(k, Hamk_bulk)
else
call ham_bulk_latticegauge (k, Hamk_bulk)
endif
endif
Hamk_bulk= (E_arc -zi* eta_arc)* ones - Hamk_bulk
call inv(Num_wann, Hamk_bulk)
do i=1, Num_wann
dos(ik, i)= aimag(Hamk_bulk(i, i))/pi
enddo
call now(time_end)
enddo
#if defined (MPI)
call mpi_allreduce(dos,dos_mpi,size(dos),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
#else
dos_mpi= dos
#endif
dos_total= eps9
do ik=1, knv3
do i=1, Num_wann
dos_total(ik)= dos_total(ik)+ dos_mpi(ik, i)
enddo
enddo
dos_selected= 0d0
do ik=1, knv3
do ig=1, NumberofSelectedOrbitals_groups
do i=1, NumberofSelectedOrbitals(ig)
dos_selected(ik, ig)= dos_selected(ik, ig)+ dos_mpi(ik, Selected_WannierOrbitals(ig)%iarray(i))
enddo
enddo
dos_unselected(ik, ig)= dos_total(ik)- dos_selected(ik, ig)
enddo
dos_max= maxval(dos_total)/4d0
outfileindex= outfileindex+ 1
if (cpuid==0)then
open(unit=outfileindex, file='fs.dat')
write(outfileindex, '(80a16, a)')'# kx', 'ky', 'kz', 'kp1', 'kp2', 'kp3', 'total A', 'A1(k,E)', 'A2(k,E)'
do ik=1, knv3
if (dos_total(ik)>dos_max) then
write(outfileindex, '(3000f16.8)')kxy_shape(:, ik)*Angstrom2atomic, kxy_plane(:, ik)*Angstrom2atomic, &
dos_total(ik), ((dos_selected(ik, ig))/dos_total(ik)*500d0+100, &
(dos_unselected(ik, ig))/dos_total(ik)*500d0+100, ig=1, NumberofSelectedOrbitals_groups)
else
write(outfileindex, '(3000f16.8)')kxy_shape(:, ik)*Angstrom2atomic, kxy_plane(:, ik)*Angstrom2atomic, &
dos_total(ik), ((dos_selected(ik, ig)), (dos_unselected(ik, ig)), &
ig=1, NumberofSelectedOrbitals_groups)
endif
if (mod(ik, nky)==0) write(outfileindex, *)' '
enddo
close(outfileindex)
endif
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
outfileindex= outfileindex+ 1
if (cpuid==0)then
open(unit=outfileindex, file='fs.dat-matlab')
write(outfileindex, '(80a16, a)')'% kx', 'ky', 'kz', 'kp1', 'kp2', 'kp3', 'total A', 'A1(k,E)', 'A2(k,E)'
write(outfileindex, '(a16, 2i16)')'% Nk1, Nk2=', Nk1, Nk2
do ik=1, knv3
write(outfileindex, '(3000f16.8)')kxy_shape(:, ik)*Angstrom2atomic, kxy_plane(:, ik)*Angstrom2atomic, &
dos_total(ik), ((dos_selected(ik, ig)), (dos_unselected(ik, ig)), &
ig=1, NumberofSelectedOrbitals_groups)
enddo
close(outfileindex)
endif
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
outfileindex= outfileindex+ 1
if (cpuid==0)then
open(unit=outfileindex, file='orbitaltexture.dat')
write(outfileindex, '(8a16, a)')'# kx', 'ky', 'kz', 'kp1', 'kp2', 'kp3', '(A1(k,E))', 'A2(k,E)'
do ik=1, knv3
if (dos_total(ik)>dos_max) then
write(outfileindex, '(300f16.8)')kxy_shape(:, ik)*Angstrom2atomic, kxy_plane(:, ik)*Angstrom2atomic, &
(dos_selected(ik, ig), dos_unselected(ik, ig), ig=1, NumberofSelectedOrbitals_groups)
endif
if (mod(ik, nky)==0) write(outfileindex, *)' '
enddo
close(outfileindex)
endif
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
zmax= maxval(log(dos_total))
zmin= minval(log(dos_total))
!> minimum and maximum value of energy bands
outfileindex= outfileindex+ 1
!> write script for gnuplot
if (cpuid==0) then
open(unit=outfileindex, file='fs.gnu')
write(outfileindex, '(a)')"set encoding iso_8859_1"
write(outfileindex, '(a)')'#set terminal postscript enhanced color'
write(outfileindex, '(a)')"#set output 'fs.eps'"
write(outfileindex, '(3a)')'set terminal pngcairo truecolor enhanced', &
' size 1920, 1680 font ",36"'
write(outfileindex, '(a)')"set output 'fs.png'"
write(outfileindex,'(a, f10.4, 2a, f10.4, a)') &
'#set palette defined ( ', zmin, ' "white", ', &
'0 "black", ', zmax,' "red" )'
write(outfileindex,'(a)') &
'set palette defined ( 0 "white", 100 "white", 101 "green", 500 "red" )'
write(outfileindex, '(a)')'#set palette rgbformulae 33,13,10'
write(outfileindex, '(a)')'unset ztics'
write(outfileindex, '(a)')'unset key'
write(outfileindex, '(a)')'set pm3d'
write(outfileindex, '(a)')'#set view equal xyz'
write(outfileindex, '(a)')'set view map'
write(outfileindex, '(a)')'set border lw 3'
write(outfileindex, '(a)')'#set xtics font ",24"'
write(outfileindex, '(a)')'#set ytics font ",24"'
write(outfileindex, '(a)')'set size ratio -1'
write(outfileindex, '(a)')'unset xtics'
write(outfileindex, '(a)')'unset ytics'
write(outfileindex, '(a)')'set colorbox'
write(outfileindex, '(a)')'set autoscale fix'
write(outfileindex, '(a)')'set pm3d interpolate 2,2'
write(outfileindex, '(2a)')"splot 'fs.dat' u 4:5:(($8)) w pm3d , \"
write(outfileindex, '(a)')" 'orbitaltexture.dat' u 4:5:(0):($7/90000):($8/90000):(0) \"
write(outfileindex, '(a)')" w vec head lw 5 lc rgb 'orange' front"
close(outfileindex)
endif
return
end subroutine orbitaltexture
subroutine fermisurface_kplane
! This subroutine calculates 3D fermi surface in the a fixed k plane
! using wannier TB method
use wmpi
use para
implicit none
integer :: ia, ik, i, j, i1, i2, ig, io
integer :: knv3, nkx, nky, ierr, nwann
real(dp) :: kz, k(3), s0(3)
! Hamiltonian of bulk system
complex(Dp), allocatable :: Hamk_bulk(:, :)
real(dp) :: time_start, time_end
real(dp) :: kxmin, kxmax, kymin, kymax
real(dp) :: zmin, zmax, dos_selected_max
real(dp) :: kxmin_shape, kxmax_shape, kymin_shape, kymax_shape
real(dp), allocatable :: kxy(:,:), kxy_shape(:,:), kxy_plane(:,:)
!> Wannier orbital projected and k-resolved spectral function
real(dp), allocatable :: dos_selected(:, :), dos_selected_mpi(:, :)
real(dp), allocatable :: dos_total(:), dos_total_mpi(:), s1(:, :)
real(dp), allocatable :: sx_selected(:, :), sy_selected(:, :), sz_selected(:, :)
real(dp), allocatable :: sx_selected_mpi(:, :), sy_selected_mpi(:, :), sz_selected_mpi(:, :)
!> spin 1/2 Pauli matrix in Wannier basis
complex(Dp),allocatable :: spin_sigma_x(:,:), spin_sigma_y(:,:), spin_sigma_z(:,:)
complex(dp), allocatable :: ones(:, :), ctemp(:, :)
logical, external :: if_in_the_WS
nkx= Nk
nky= Nk
knv3= nkx*nky
allocate( kxy(3, knv3))
allocate( kxy_shape(3, knv3))
allocate( kxy_plane(3, knv3))
kxy=0d0
kxy_shape=0d0
kxy_plane=0d0
ik =0
do i= 1, nkx
do j= 1, nky
ik =ik +1
kxy(:, ik)= K3D_start+ K3D_vec1*(i-1)/dble(nkx-1)+ K3D_vec2*(j-1)/dble(nky-1) &
-(K3D_vec1+K3D_vec2)/2d0
kxy_shape(:, ik)= kxy(1, ik)* Origin_cell%Kua+ kxy(2, ik)* Origin_cell%Kub+ kxy(3, ik)* Origin_cell%Kuc
call rotate_k3_to_kplane(kxy_shape(:, ik), kxy_plane(:, ik))
enddo
enddo
i1=1
i2=2
kxmin_shape=minval(kxy_shape(i1,:))
kxmax_shape=maxval(kxy_shape(i1,:))
kymin_shape=minval(kxy_shape(i2,:))
kymax_shape=maxval(kxy_shape(i2,:))
allocate(Hamk_bulk(Num_wann, Num_wann))
allocate( dos_total(knv3), dos_total_mpi(knv3))
allocate( dos_selected (knv3, NumberofSelectedOrbitals_groups))
allocate( dos_selected_mpi(knv3, NumberofSelectedOrbitals_groups))
dos_total= eps12; dos_total_mpi= eps12
dos_selected= eps12; dos_selected_mpi= eps12
if (SOC>0 .and. BulkSpintexture_calc) then
allocate( sx_selected(knv3, NumberofSelectedOrbitals_groups))
allocate( sy_selected(knv3, NumberofSelectedOrbitals_groups))
allocate( sz_selected(knv3, NumberofSelectedOrbitals_groups))
allocate( sx_selected_mpi(knv3, NumberofSelectedOrbitals_groups))
allocate( sy_selected_mpi(knv3, NumberofSelectedOrbitals_groups))
allocate( sz_selected_mpi(knv3, NumberofSelectedOrbitals_groups))
allocate(spin_sigma_x(Num_wann,Num_wann), spin_sigma_y(Num_wann,Num_wann), spin_sigma_z(Num_wann,Num_wann))
sx_selected= 0d0; sx_selected_mpi= 0d0
sy_selected= 0d0; sy_selected_mpi= 0d0
sz_selected= 0d0; sz_selected_mpi= 0d0
spin_sigma_x= 0d0; spin_sigma_y= 0d0; spin_sigma_z= 0d0
endif
if (SOC>0 .and. BulkSpintexture_calc) then
if (index(Particle,'phonon')/=0) then
stop "ERROR: we don't support spintexture calculation for phonon system"
endif
Nwann= Num_wann/2
!> spin operator matrix
!> this part is package dependent.
!if (index( Package, 'VASP')/=0.or. index( Package, 'Wien2k')/=0 &
! .or. index( Package, 'Abinit')/=0.or. index( Package, 'openmx')/=0) then
do j=1, Nwann
spin_sigma_x(j, Nwann+j)=1.0d0
spin_sigma_x(j+Nwann, j)=1.0d0
spin_sigma_y(j, Nwann+j)=-zi
spin_sigma_y(j+Nwann, j)=zi
spin_sigma_z(j, j)= 1d0
spin_sigma_z(j+Nwann, j+Nwann)=-1d0
enddo
!else
! if (cpuid.eq.0) write(stdout, *)'Error: please report your software generating tight binding and wannier90.wout to me'
! if (cpuid.eq.0) write(stdout, *)'wuquansheng@gmail.com'
! stop 'Error: please report your software and wannier90.wout to wuquansheng@gmail.com'
!endif
endif
allocate(ones(Num_wann, Num_wann), ctemp(Num_wann, Num_wann))
ones= 0d0
do i=1, Num_wann
ones(i, i)= 1d0
enddo
time_start= 0d0
time_end= 0d0
do ik= 1+cpuid, knv3, num_cpu
if (cpuid==0.and. mod(ik/num_cpu, 500)==0) &
write(stdout, *) 'FS_Plane, ik ', ik, 'knv3',knv3, 'time left', &
(knv3-ik)*(time_end- time_start)/num_cpu, ' s'
call now(time_start)
k = kxy(:, ik)
!> if the k points is not in the first Wigner-Seitz cell, then ignore it
if (Translate_to_WS_calc)then
if (.not.if_in_the_WS(k)) then
cycle
endif
endif
! calculation bulk hamiltonian
Hamk_bulk= 0d0
call ham_bulk_latticegauge(k, Hamk_bulk)
Hamk_bulk= (E_arc -zi* eta_arc)* ones - Hamk_bulk
!> Hamk_bulk is the Green's function after subroutine inv.
call inv(Num_wann, Hamk_bulk)
do i=1, Num_wann
dos_total(ik)= dos_total(ik)+ aimag(Hamk_bulk(i, i))/pi
enddo
do ig=1, NumberofSelectedOrbitals_groups
do i=1, NumberofSelectedOrbitals(ig)
io=Selected_WannierOrbitals(ig)%iarray(i)
dos_selected(ik, ig)= dos_selected(ik, ig)+ aimag(Hamk_bulk(io, io))/pi
enddo
enddo
call now(time_end)
if (SOC>0.and.BulkSpintexture_calc) then
!>> calculate spin-resolved bulk energy spectrum
!> spin x
call mat_mul(Num_wann, Hamk_bulk, spin_sigma_x, ctemp)
do ig=1, NumberofSelectedOrbitals_groups
do i=1, NumberofSelectedOrbitals(ig)
io=Selected_WannierOrbitals(ig)%iarray(i)
sx_selected(ik, ig)= sx_selected(ik, ig)+ aimag(ctemp(io, io))/pi
enddo
enddo
!> spin y
call mat_mul(Num_wann, Hamk_bulk, spin_sigma_y, ctemp)
do ig=1, NumberofSelectedOrbitals_groups
do i=1, NumberofSelectedOrbitals(ig)
io=Selected_WannierOrbitals(ig)%iarray(i)
sy_selected(ik, ig)= sy_selected(ik, ig)+ aimag(ctemp(io, io))/pi
enddo
enddo
!> spin z
call mat_mul(Num_wann, Hamk_bulk, spin_sigma_z, ctemp)
do ig=1, NumberofSelectedOrbitals_groups
do i=1, NumberofSelectedOrbitals(ig)
io=Selected_WannierOrbitals(ig)%iarray(i)
sz_selected(ik, ig)= sz_selected(ik, ig)+ aimag(ctemp(io, io))/pi
enddo
enddo
endif
enddo
#if defined (MPI)
call mpi_allreduce(dos_total,dos_total_mpi,size(dos_total),&
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)
if (SOC>0.and.BulkSpintexture_calc) then
call mpi_allreduce(sx_selected,sx_selected_mpi,size(sx_selected),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
call mpi_allreduce(sy_selected,sy_selected_mpi,size(sy_selected),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
call mpi_allreduce(sz_selected,sz_selected_mpi,size(sz_selected),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
endif
#else
dos_total_mpi= dos_total
dos_selected_mpi= dos_selected
if (SOC>0.and.BulkSpintexture_calc) then
sx_selected_mpi= sx_selected
sy_selected_mpi= sy_selected
sz_selected_mpi= sz_selected
endif
#endif
if (SOC>0.and.BulkSpintexture_calc) then
do ig=1, NumberofSelectedOrbitals_groups
dos_selected_max= maxval(dos_selected_mpi(:, ig))
do ik=1, knv3
if (dos_selected_mpi(ik, ig)<dos_selected_max/20d0) then
sx_selected_mpi(ik, ig)= eps12
sy_selected_mpi(ik, ig)= eps12
sz_selected_mpi(ik, ig)= eps12
endif
enddo
enddo
allocate(s1(3, NumberofSelectedOrbitals_groups))
s1= 0d0
endif
outfileindex= outfileindex+ 1
if (SOC>0.and.BulkSpintexture_calc) then
if (cpuid==0)then
open(unit=outfileindex, file='bulkspintext.dat')
write(outfileindex, &
"('#', a12, 5a16, 3X, '| A(k,E)', a6, 100(4X,'|group ', i2, ': A ', 13X, 'sx', 14X, 'sy', 14X, 'sz'))")&
'kx', 'ky', 'kz', 'kp1', 'kp2', 'kp3', 'total',&
(i, i=1, NumberofSelectedOrbitals_groups)
write(outfileindex, "('#column', i5, 3000i16)")(i, i=1, 7+4*NumberofSelectedOrbitals_groups)
do ik=1, knv3
do ig=1, NumberofSelectedOrbitals_groups
s0(1)= sx_selected_mpi(ik, ig) /dos_selected_mpi(ik, ig)
s0(2)= sy_selected_mpi(ik, ig) /dos_selected_mpi(ik, ig)
s0(3)= sz_selected_mpi(ik, ig) /dos_selected_mpi(ik, ig)
call rotate_k3_to_kplane(s0, s1(:, ig))
enddo
write(outfileindex, '(3000E16.8)')kxy_shape(:, ik)*Angstrom2atomic, kxy_plane(:, ik)*Angstrom2atomic, &
dos_total_mpi(ik), (dos_selected_mpi(ik, ig), s1(:, ig), ig=1, NumberofSelectedOrbitals_groups)
if (mod(ik, nky)==0) write(outfileindex, *)' '
enddo
close(outfileindex)
endif
else
if (cpuid==0)then
open(unit=outfileindex, file='fs_kplane.dat')
write(outfileindex, "('#', a12, 5a16, 3X, '| A(k,E)', a6, 100(8X,'group ', i2))")&
'kx', 'ky', 'kz', 'kp1', 'kp2', 'kp3', 'total',&
(i, i=1, NumberofSelectedOrbitals_groups)
write(outfileindex, "('#column', i5, 3000i16)")(i, i=1, 7+NumberofSelectedOrbitals_groups)
do ik=1, knv3
write(outfileindex, '(3000E16.8)')kxy_shape(:, ik)*Angstrom2atomic, kxy_plane(:, ik)*Angstrom2atomic, &
dos_total_mpi(ik), (dos_selected_mpi(ik, ig), ig=1, NumberofSelectedOrbitals_groups)
if (mod(ik, nky)==0) write(outfileindex, *)' '
enddo
close(outfileindex)
endif
endif
#if defined (MPI)
call mpi_barrier(mpi_cmw, ierr)
#endif
!> minimum and maximum value of dos_total
zmax= maxval(log(dos_total_mpi))
zmin= minval(log(dos_total_mpi))
outfileindex= outfileindex+ 1
!> write script for gnuplot
if (cpuid==0) then
if (SOC>0.and.BulkSpintexture_calc) then
open(unit=outfileindex, file='bulkspintext.gnu')
write(outfileindex, '(a)')"set encoding iso_8859_1"
write(outfileindex, '(3a)')'set terminal pngcairo truecolor enhanced', &
' size 1920, 1680 font ",36"'
write(outfileindex, '(a)')"set output 'bulkspintext.png'"
write(outfileindex,'(a, f10.4, a, f10.4, a, f10.4, a)') &
'set palette defined ( ', zmin, ' "white", ', &
(zmin+zmax)/2d0,' "black", ', zmax,' "red" )'
write(outfileindex, '(a)')'#set palette rgbformulae 33,13,10'
write(outfileindex, '(a)')'unset ztics'
write(outfileindex, '(a)')'unset key'
write(outfileindex, '(a)')'set pm3d'
write(outfileindex, '(a)')'#set view equal xyz'
write(outfileindex, '(a)')'set view map'
write(outfileindex, '(a)')'set border lw 3'
write(outfileindex, '(a)')'#set xtics font ",24"'
write(outfileindex, '(a)')'#set ytics font ",24"'
write(outfileindex, '(a)')'set size ratio -1'
write(outfileindex, '(a)')'unset xtics'
write(outfileindex, '(a)')'unset ytics'
write(outfileindex, '(a)')'set colorbox'
write(outfileindex, '(a)')'set autoscale fix'
write(outfileindex, '(a)')'set pm3d interpolate 2,2'
write(outfileindex, '(2a)')"splot 'bulkspintext.dat' u 4:5:(log($8+0.1)) w pm3d, \"
write(outfileindex, '(a)')" 'bulkspintext.dat' u 4:5:(0):($9/5.00):($10/5.00):(0) w vec head lw 5 lc rgb 'orange' front"
close(outfileindex)
else
open(unit=outfileindex, file='fs_kplane.gnu')
write(outfileindex, '(a)')"set encoding iso_8859_1"
write(outfileindex, '(3a)')'set terminal pngcairo truecolor enhanced', &
' size 1920, 1680 font ",36"'
write(outfileindex, '(a)')"set output 'fs_kplane.png'"
write(outfileindex,'(a, f10.4, a, f10.4, a, f10.4, a)') &
'set palette defined ( ', zmin, ' "white", ', &
(zmin+zmax)/2d0,' "black", ', zmax,' "red" )'
write(outfileindex, '(a)')'#set palette rgbformulae 33,13,10'
write(outfileindex, '(a)')'unset ztics'
write(outfileindex, '(a)')'unset key'
write(outfileindex, '(a)')'set pm3d'
write(outfileindex, '(a)')'#set view equal xyz'
write(outfileindex, '(a)')'set view map'
write(outfileindex, '(a)')'set border lw 3'
write(outfileindex, '(a)')'#set xtics font ",24"'
write(outfileindex, '(a)')'#set ytics font ",24"'
write(outfileindex, '(a)')'set size ratio -1'
write(outfileindex, '(a)')'unset xtics'
write(outfileindex, '(a)')'unset ytics'
write(outfileindex, '(a)')'set colorbox'
write(outfileindex, '(a)')'set autoscale fix'
write(outfileindex, '(a)')'set pm3d interpolate 2,2'
write(outfileindex, '(2a)')"splot 'fs_kplane.dat' u 4:5:(log($7+0.1)) w pm3d"
close(outfileindex)
endif
endif
return
end subroutine fermisurface_kplane
subroutine fermisurface_stack
!> This subroutine calculates 3D fermi surfaces in the a fixed k plane
!> using wannier TB method
!> With the given k-plane set by KPLANE_BULK, we calculate a set of spectrum with the kplane
!> parallel to the KPLANE_BULK
use wmpi
use para
implicit none
integer :: ia, ik, i, j, i1, i2, ig
integer :: knv3, nkx, nky, ierr, ikz
real(dp) :: kz(3), k(3), k1(3), k2(3)
! Hamiltonian of bulk system
complex(Dp), allocatable :: Hamk_bulk(:, :)
real(dp) :: time_start, time_end
real(dp) :: kxmin, kxmax, kymin, kymax
real(dp) :: zmin, zmax
real(dp) :: kxmin_shape, kxmax_shape, kymin_shape, kymax_shape
real(dp), allocatable :: kxy(:,:), kxy_shape(:,:), kxy_plane(:,:)
!> atomic and k resolved spectral function
real(dp), allocatable :: dos(:, :), dos_mpi(:, :), dos_atom(:, :)
real(dp), allocatable :: dos_atom_mpi(:, :), dos_total(:, :)
complex(dp), allocatable :: ones(:,:)
integer, allocatable :: orbitals_end(:)
logical, external :: if_in_the_WS
real(dp), external :: norm
allocate(Hamk_bulk(Num_wann, Num_wann))
allocate(orbitals_end(Origin_cell%Num_atoms))
orbitals_end(1)= Origin_cell%nprojs(1)
do ia=2, Origin_cell%Num_atoms
orbitals_end(ia)= orbitals_start(ia)+ Origin_cell%nprojs(ia)- 1
enddo
nkx= Nk
nky= Nk
allocate( kxy(3, nkx*nky))
allocate( kxy_shape(3, nkx*nky))
allocate( kxy_plane(3, nkx*nky))
kxy=0d0
kxy_shape=0d0
kxy_plane=0d0
!> a k-plane defined by KPLANE_BULK
ik =0
do i= 1, nkx
do j= 1, nky
ik =ik +1
kxy(:, ik)= K3D_start+ K3D_vec1*(i-1)/dble(nkx-1)+ K3D_vec2*(j-1)/dble(nky-1) &
-(K3D_vec1+K3D_vec2)/2d0
kxy_shape(:, ik)= kxy(1, ik)* Origin_cell%Kua+ kxy(2, ik)* Origin_cell%Kub+ kxy(3, ik)* Origin_cell%Kuc
call rotate_k3_to_kplane(kxy_shape(:, ik), kxy_plane(:, ik))
enddo
enddo
!> find a vector that perpendicular to the k-plane
call direct_cart_rec(K3D_vec1, k1)
call direct_cart_rec(K3D_vec2, k2)
call cross_product(k1, k2, kz)
kz=kz/norm(kz)*2d0
k2=kz
call cart_direct_rec(k2, kz)
i1=1
i2=2
kxmin_shape=minval(kxy_shape(i1,:))
kxmax_shape=maxval(kxy_shape(i1,:))
kymin_shape=minval(kxy_shape(i2,:))
kymax_shape=maxval(kxy_shape(i2,:))
knv3= nkx*nky
allocate( dos_atom (knv3, Origin_cell%Num_atoms))
allocate( dos_atom_mpi(knv3, Origin_cell%Num_atoms))
allocate( dos (knv3, Num_wann))
allocate( dos_mpi(knv3, Num_wann))
allocate( dos_total(knv3, NumberofSelectedOrbitals_groups))
dos = 0d0
dos_mpi= 0d0
dos_atom = 0d0
dos_atom_mpi= 0d0
allocate(ones(Num_wann, Num_wann))
ones= 0d0
do i=1, Num_wann
ones(i, i)= 1d0
enddo
time_start= 0d0
time_end= 0d0
do ik= 1+cpuid, knv3, num_cpu
if (cpuid==0.and. mod(ik/num_cpu, 500)==0) &
write(stdout, *) 'FS_Plane, ik ', ik, 'knv3',knv3, 'time left', &
(knv3-ik)*(time_end- time_start)/num_cpu, ' s'
call now(time_start)
k = kxy(:, ik)
if (Translate_to_WS_calc)then
if (.not.if_in_the_WS(k)) then
dos(ik, :)= eps6
dos_atom(ik, :)= eps6
cycle
endif
endif
do ikz=1, Nk3
k = kxy(:, ik)+ kz*(ikz-1)/Nk3
! calculation bulk hamiltonian
Hamk_bulk= 0d0
call ham_bulk_latticegauge(k, Hamk_bulk)
Hamk_bulk= (E_arc -zi* eta_arc)* ones - Hamk_bulk
call inv(Num_wann, Hamk_bulk)
do i=1, Num_wann
dos(ik, i)=dos(ik, i)+ aimag(Hamk_bulk(i, i))/pi
enddo
do ia=1, Origin_cell%Num_atoms
do i=orbitals_start(ia), orbitals_end(ia)
dos_atom(ik, ia)= dos_atom(ik, ia)+ aimag(Hamk_bulk(i, i))/pi
enddo
if (SOC>0) then
do i=orbitals_start(ia)+num_wann/2, orbitals_end(ia)+num_wann/2
dos_atom(ik, ia)= dos_atom(ik, ia)+ aimag(Hamk_bulk(i, i))/pi
enddo
endif
enddo ! ia
enddo ! ikz
call now(time_end)
enddo
#if defined (MPI)
call mpi_allreduce(dos,dos_mpi,size(dos),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
call mpi_allreduce(dos_atom,dos_atom_mpi,size(dos_atom),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
#else
dos_mpi= dos
dos_atom_mpi= dos_atom
#endif
dos_total= eps9
do ik=1, knv3
do ig=1, NumberofSelectedOrbitals_groups
do i=1, NumberofSelectedOrbitals(ig)
dos_total(ik, ig)= dos_total(ik, ig)+ dos_mpi(ik, Selected_WannierOrbitals(ig)%iarray(i))
enddo
enddo
enddo
outfileindex= outfileindex+ 1
if (cpuid==0)then
open(unit=outfileindex, file='fs_stack.dat')
write(outfileindex, '(7a16, a)')'# kx', 'ky', 'kz', 'kp1', 'kp2', 'kp3', '(A(k,E))', 'A(i,k,E),i=1,numatoms'
do ik=1, knv3
write(outfileindex, '(3000f16.8)')kxy_shape(:, ik)*Angstrom2atomic, kxy_plane(:, ik)*Angstrom2atomic, &
(dos_total(ik, ig), ig=1, NumberofSelectedOrbitals_groups), (dos_atom_mpi(ik, :))
if (mod(ik, nky)==0) write(outfileindex, *)' '
enddo
close(outfileindex)
endif
zmax= maxval(log(dos_total))
zmin= minval(log(dos_total))
!> minimum and maximum value of energy bands
outfileindex= outfileindex+ 1
!> write script for gnuplot
if (cpuid==0) then
open(unit=outfileindex, file='fs_stack.gnu')
write(outfileindex, '(a)')"set encoding iso_8859_1"
write(outfileindex, '(a)')'#set terminal postscript enhanced color'
write(outfileindex, '(a)')"#set output 'fs_stack.eps'"
write(outfileindex, '(3a)')'set terminal pngcairo truecolor enhanced', &
' size 1920, 1680 font ",36"'
write(outfileindex, '(a)')"set output 'fs_stack.png'"
write(outfileindex,'(a, f10.4, 2a, f10.4, a)') &
'set palette defined ( ', zmin, ' "white", ', &
'0 "black", ', zmax,' "red" )'
write(outfileindex, '(a)')'#set palette rgbformulae 33,13,10'
write(outfileindex, '(a)')'unset ztics'
write(outfileindex, '(a)')'unset key'
write(outfileindex, '(a)')'set pm3d'
write(outfileindex, '(a)')'#set view equal xyz'
write(outfileindex, '(a)')'set view map'
write(outfileindex, '(a)')'set border lw 3'
write(outfileindex, '(a)')'#set xtics font ",24"'
write(outfileindex, '(a)')'#set ytics font ",24"'
write(outfileindex, '(a)')'set size ratio -1'
write(outfileindex, '(a)')'unset xtics'
write(outfileindex, '(a)')'unset ytics'
write(outfileindex, '(a)')'set colorbox'
!write(outfileindex, '(a, f10.5, a, f10.5, a)')'set xrange [', kxmin, ':', kxmax, ']'
!write(outfileindex, '(a, f10.5, a, f10.5, a)')'set yrange [', kymin, ':', kymax, ']'
write(outfileindex, '(a)')'set autoscale fix'
!write(outfileindex, '(a, f10.5, a, f10.5, a)')'set xrange [', kxmin_shape, ':', kxmax_shape, ']'
!write(outfileindex, '(a, f10.5, a, f10.5, a)')'set yrange [', kymin_shape, ':', kymax_shape, ']'
write(outfileindex, '(a)')'set pm3d interpolate 2,2'
write(outfileindex, '(2a)')"splot 'fs_stack.dat' u 4:5:(log($7+0.1)) w pm3d"
close(outfileindex)
endif
return
end subroutine fermisurface_stack
subroutine gapshape3D
! This subroutine get the k points in the BZ at which the gap is smaller then
! Gap_threshold
use wmpi
use para
implicit none
integer :: ik, i, j, l
integer :: knv3