-
Notifications
You must be signed in to change notification settings - Fork 0
/
ek_bulk_polar.f90
executable file
·315 lines (269 loc) · 8.93 KB
/
ek_bulk_polar.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
subroutine ek_bulk_polar
! Calculate bulk's energy band using wannier TB method.
! This subroutine is very useful for nodal-line materials.
! You have to define the center, k plane and radius
! Copyright (c) 2010 QuanSheng Wu. All rights reserved.
use wmpi
use para
implicit none
integer :: ik, i, j, ia, knv3
integer :: ierr, nktheta, nkr, nwan
real(Dp) :: k (3), k1(3)
real(dp) :: kz, ktheta, kr, krmax, mingap
real(Dp) :: W(Num_wann)
real(dp), allocatable :: k12(:,:), k123(:,:), krtheta(:,:)
!> center
real(dp) :: center_direct(3)
real(dp) :: center_cart(3)
! Hamiltonian of bulk system
complex(Dp) :: Hamk_bulk(Num_wann,Num_wann)
! eigen value of H
real(dp), allocatable :: gap(:,:)
real(dp), allocatable :: eigv(:,:)
real(dp), allocatable :: eigv_mpi(:,:)
!> pauli matrix in basis
complex(dp), allocatable :: sigmax(:, :)
complex(dp), allocatable :: sigmay(:, :)
complex(dp), allocatable :: sigmaz(:, :)
complex(dp), allocatable :: mat1(:, :)
complex(dp), allocatable :: eigvec_conj(:, :)
!> spin expect value for each band and each k point
real(dp), allocatable :: sx(:, :)
real(dp), allocatable :: sy(:, :)
real(dp), allocatable :: sz(:, :)
nktheta= 100
nkr= Nk
knv3= nktheta*nkr
allocate( k12(2, knv3))
allocate( k123(3, knv3))
allocate( krtheta(2, knv3))
allocate( gap (4, nktheta))
allocate( eigv (4, knv3))
allocate( eigv_mpi(4, knv3))
allocate( eigvec_conj(Num_wann, Num_wann))
allocate( sigmax(Num_wann, Num_wann))
allocate( sigmay(Num_wann, Num_wann))
allocate( sigmaz(Num_wann, Num_wann))
allocate( mat1(Num_wann, Num_wann))
allocate( sx(Num_wann, knv3))
allocate( sy(Num_wann, knv3))
allocate( sz(Num_wann, knv3))
k12 = 0d0
k123= 0d0
krtheta = 0d0
eigv = 0d0
eigv_mpi= 0d0
krmax= 0.6d0
kz= 0d0
ik =0
!> for IrF DFT
!> xz plane, center is (0.0, 0.5, 0.5)
!> for IrF TB
!> xz plane, center is (0.5, 0.0, 0.5)
!center_direct= (/ 0.5d0, 0.0d0, 0.5d0/)
!call direct_cart_rec(center_direct, center_cart)
!ik = 0
!do i= 1, nktheta
! ktheta= 2d0*pi*(i-1)/dble(nktheta)
! do j= 1, nkr
! kr= krmax*(j-1)/dble(nkr)
! ik =ik +1
! k12(1, ik)= kr*cos(ktheta)+ center_cart(1)
! k12(2, ik)= kr*sin(ktheta)+ center_cart(3)
! krtheta(1, ik)= ktheta
! krtheta(2, ik)= kr
! k123(1, ik) = k12(1, ik)
! k123(2, ik)= center_cart(2)
! k123(3, ik) = k12(2, ik)
! enddo
!enddo
!> For DFT yz plane, center is (0.5, 0.0, 0.5)
!center_direct= (/ 0.5d0, 0.0d0, 0.5d0/)
!call direct_cart_rec(center_direct, center_cart)
!ik = 0
!do i= 1, nktheta
! ktheta= 2d0*pi*(i-1)/dble(nktheta)
! do j= 1, nkr
! kr= krmax*(j-1)/dble(nkr)
! ik =ik +1
! k12(1, ik)= kr*cos(ktheta)+ center_cart(2)
! k12(2, ik)= kr*sin(ktheta)+ center_cart(3)
! krtheta(1, ik)= ktheta
! krtheta(2, ik)= kr
! k123(1, ik)= center_cart(1)
! k123(2, ik) = k12(1, ik)
! k123(3, ik) = k12(2, ik)
! enddo
!enddo
!>> for TB
!> xy plane, center is (0.0, 0.0, 0.0)
!center_direct= (/ 0.0d0, 0.0d0, -0.0d0/)
!call direct_cart_rec(center_direct, center_cart)
!ik = 0
!do i= 1, nktheta
! ktheta= 2d0*pi*(i-1)/dble(nktheta)
! do j= 1, nkr
! kr= krmax*(j-1)/dble(nkr)
! ik =ik +1
! k12(1, ik)= kr*cos(ktheta)+ center_cart(1)
! k12(2, ik)= kr*sin(ktheta)+ center_cart(2)
! krtheta(1, ik)= ktheta
! krtheta(2, ik)= kr
! k123(1, ik) = k12(1, ik)
! k123(2, ik) = k12(2, ik)
! k123(3, ik)= center_cart(3)
! enddo
!enddo
!> xz plane, center is (0.5, 0.0, 0.5)
!center_direct= (/ 0.5d0, 0.0d0, 0.5d0/)
!center_direct= (/ 0.5d0, 0.0d0, -0.5d0/)
!call direct_cart_rec(center_direct, center_cart)
!ik = 0
!do i= 1, nktheta
! ktheta= 2d0*pi*(i-1)/dble(nktheta)
! do j= 1, nkr
! kr= krmax*(j-1)/dble(nkr)
! ik =ik +1
! k12(1, ik)= kr*cos(ktheta)+ center_cart(1)
! k12(2, ik)= kr*sin(ktheta)+ center_cart(3)
! krtheta(1, ik)= ktheta
! krtheta(2, ik)= kr
! k123(1, ik) = k12(1, ik)
! k123(2, ik)= center_cart(2)
! k123(3, ik) = k12(2, ik)
! enddo
!enddo
!> For yz plane, center is (0.0, 0.5, 0.5)
!center_direct= (/ 0.0d0,-0.5d0,-0.5d0/)
!center_direct= (/ 1.0d0, 0.5d0, 0.5d0/)
!call direct_cart_rec(center_direct, center_cart)
!ik = 0
!do i= 1, nktheta
! ktheta= 2d0*pi*(i-1)/dble(nktheta)
! do j= 1, nkr
! kr= krmax*(j-1)/dble(nkr)
! ik =ik +1
! k12(1, ik)= kr*cos(ktheta)+ center_cart(2)
! k12(2, ik)= kr*sin(ktheta)+ center_cart(3)
! krtheta(1, ik)= ktheta
! krtheta(2, ik)= kr
! k123(1, ik)= center_cart(1)
! k123(2, ik) = k12(1, ik)
! k123(3, ik) = k12(2, ik)
! enddo
!enddo
!>> for CuTlTe2
!> 110 plane, center is X (0.0, 0.0, 0.5)
center_direct= (/ 0.0d0, 0.0d0, 0.5d0/)
call direct_cart_rec(center_direct, center_cart)
if (cpuid.eq.0) write(*, '(10f8.5)') center_direct
if (cpuid.eq.0) write(*, '(10f8.5)') center_cart
ik = 0
krmax= 0.3d0
do i= 1, nktheta
ktheta= 2d0*pi*(i-1)/dble(nktheta)
do j= 1, nkr
kr= krmax*(j-1)/dble(nkr)
ik =ik +1
k123(1, ik)= kr*cos(ktheta)+ center_cart(1)
k123(2, ik)= kr*cos(ktheta)+ center_cart(2)
k123(3, ik)= kr*sin(ktheta)+ center_cart(3)
enddo
enddo
!>> for CuTlTe2
!> 1-10 plane, center is U (0.5, 0.5, 0.0)
!center_direct= (/ 0.5d0, 0.5d0, 0.0d0/)
!call direct_cart_rec(center_direct, center_cart)
!if (cpuid.eq.0) write(*, '(10f8.5)') center_direct
!if (cpuid.eq.0) write(*, '(10f8.5)') center_cart
!ik = 0
!krmax= 0.3d0
!do i= 1, nktheta
! ktheta= 2d0*pi*(i-1)/dble(nktheta)
! do j= 1, nkr
! kr= krmax*(j-1)/dble(nkr)
! ik =ik +1
! k123(1, ik)= kr*cos(ktheta)+ center_cart(1)
! k123(2, ik)= -kr*cos(ktheta)+ center_cart(2)
! k123(3, ik)= kr*sin(ktheta)+ center_cart(3)
! enddo
!enddo
!> set Pauli matrix
sigmax= 0d0
sigmay= 0d0
sigmaz= 0d0
if (soc>0) then
nwan= Num_wann/2
else
stop 'calculate spin texture need soc'
endif
i= 0
do ia=1, Origin_cell%Num_atoms
do j=1, Origin_cell%nprojs(ia)
i= i+ 1
sigmax(i, i+nwan) = 1d0
sigmax(i+nwan, i) = 1d0
sigmay(i, i+nwan) = -zi
sigmay(i+nwan, i) = zi
sigmaz(i, i) = 1d0
sigmaz(i+nwan, i+nwan) = -1d0
enddo !> number of projectors
enddo ! ia number of atoms
if (Numoccupied> Num_wann ) stop ' Numoccupied should less than Num_wann'
do ik= 1+cpuid, knv3, num_cpu
if (cpuid==0) print * , 'ik' , ik
k1= k123(:, ik)
!> from cartisen coordinate to direct coordinate
call cart_direct_rec(k1, k)
! calculation bulk hamiltonian
Hamk_bulk= 0d0
call ham_bulk_atomicgauge(k, Hamk_bulk)
!> diagonalization by call zheev in lapack
W= 0d0
call eigensystem_c( 'V', 'U', Num_wann ,Hamk_bulk, W)
eigvec_conj= conjg(transpose(Hamk_bulk))
eigv(:, ik)= W(Numoccupied-1:Numoccupied+2)
call mat_mul(Num_wann,eigvec_conj, sigmax, mat1)
call mat_mul(Num_wann,mat1, Hamk_bulk, sigmax)
call mat_mul(Num_wann,eigvec_conj, sigmay, mat1)
call mat_mul(Num_wann,mat1, Hamk_bulk, sigmay)
call mat_mul(Num_wann,eigvec_conj, sigmaz, mat1)
call mat_mul(Num_wann,mat1, Hamk_bulk, sigmaz)
do i=1, Num_wann
sx(i, ik)= sigmax(i, i)
sy(i, ik)= sigmay(i, i)
sz(i, ik)= sigmaz(i, i)
enddo
enddo ! ik
#if defined (MPI)
call mpi_allreduce(eigv,eigv_mpi,size(eigv),&
mpi_dp,mpi_sum,mpi_cmw,ierr)
#else
eigv_mpi= eigv
#endif
eigv= eigv_mpi
ik= 0
do i=1, nktheta
mingap= 999d0
do j= 1, nkr
ik= ik+ 1
if ((eigv(3, ik)- eigv(2, ik))< mingap) then
mingap= eigv(3, ik)- eigv(2, ik)
gap(1, i)= mingap
gap(2:4, i)= k123(:, ik)
endif
enddo
enddo
outfileindex= outfileindex+ 1
if (cpuid==0)then
open(unit=outfileindex, file='bulkek_polar.dat')
do i=1, nktheta
write(outfileindex, '(1000f19.9)')gap(:, i)
enddo
do i=1, knv3
! write(outfileindex, '(1000f19.9)')k123(:, i)
enddo
close(outfileindex)
endif
return
end subroutine ek_bulk_polar