/
reduk.f90
324 lines (292 loc) · 10 KB
/
reduk.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
! The code was developed at the Fritz Haber Institute, and
! the intellectual properties and copyright of this file
! are with the Max Planck Society. When you use it, please
! cite R. Gomez-Abal, X. Li, C. Ambrosch-Draxl, M. Scheffler,
! Extended linear tetrahedron method for the calculation of q-dependent
! dynamical response functions, to be published in Comp. Phys. Commun. (2010)
!BOP
!
! !ROUTINE: reduk
!
! !INTERFACE:
subroutine reduk(nsymt,divsh,weight)
!
! !DESCRIPTION:
!
! This subroutine calculates a set of reduced k-points for the integration
! with the tetrahedron method. The k-points are reduced by the symmetry and
! a weight is put to each irreducible k-point to tell how many points it
! represents before the symmetry operation.
!
! !USES:
use kgen_internals
!<sag>
use control, only: tetraifc, tetradbglv
!</sag>
implicit none
!<sag>
real(8), parameter :: epslat=1.d-6
real(8) :: kpr(3),kprt(3)
logical, allocatable :: done(:)
!</sag>
! !INPUT PARAMETERS:
integer(4), intent(in) :: nsymt ! Number of symmetry operations
integer(4), intent(in) :: divsh
! !OUTPUT PARAMETERS:
integer(4), intent(out):: weight(*) ! weight of the
! reducible k-point
! !LOCAL VARIABLES:
integer(4), allocatable :: kpav(:)
integer(4) :: i,j,i1,i2,i3
integer(4) :: onek(3),nkp(3),kk(3),kpid,nkpid,ktp(3)
integer(4) :: members,starr(48)
integer(4), external :: jget
integer(4), external :: idkp
!EOP
!BOC
external jset
allocate(kpav(div(1)*div(2)*div(3)/32+1))
nirkp=0
do i=1,div(1)*div(2)*div(3)
call jset(kpav,i,0)
enddo
! write(6,*)'------------------------------------------------------'
! write(6,*)' reduk: begin'
! write(6,*)'------------------------------------------------------'
!<sag>
if (trim(tetraifc)=='wien2k') then
!!! *** OLD CODE WITH BUGS ***
!!$ ! original code
!!$ do i1=0,div(1)-1
!!$ kk(1)=divsh*i1+shift(1)
!!$ onek(1)=i1
!!$ do i2=0,div(2)-1
!!$ kk(2)=divsh*i2+shift(2)
!!$ onek(2)=i2
!!$ do i3=0,div(3)-1
!!$ kk(3)=divsh*i3+shift(3)
!!$ onek(3)=i3
!!$ kpid=idkp(onek)
!!$ if(jget(kpav,kpid).eq.0) then
!!$ nirkp=nirkp+1
!!$ ikpid(kpid)=nirkp
!!$ do i=1,48
!!$ starr(i)=0
!!$ enddo
!!$ ! write(24,*)
!!$ ! write(24,'(i6,3i4)')divsh,div
!!$ ! write(24,'(3i4," ",3i4)')onek,kk
!!$ do i=1,nsymt
!!$ do j=1,3
!!$ nkp(j)=mod(iio(j,1,i)*kk(1)+iio(j,2,i)*kk(2)+&
!!$ & iio(j,3,i)*kk(3),divsh*div(j))
!!$ ktp(j)=nkp(j)
!!$ nkp(j)=nkp(j)+(1-isign(1,nkp(j)))*divsh*div(j)/2
!!$ nkp(j)=(nkp(j)-shift(j))/divsh
!!$ enddo
!!$ ! write(24,'(i6,3i4," ",3i4)')i,ktp,nkp
!!$ nkpid=idkp(nkp)
!!$ ! write(24,'(4x,2i4)')nirkp,nkpid
!!$ call jset(kpav,nkpid,1)
!!$ redkp(nkpid)=kpid
!!$ starr(i)=nkpid
!!$ enddo
!!$ members=0
!!$ do i=1,nsymt
!!$ if(starr(i).ne.0) then
!!$ members=members+1
!!$ do j=i+1,nsymt
!!$ if(starr(i).eq.starr(j)) starr(j)=0
!!$ enddo
!!$ endif
!!$ enddo
!!$ weight(nirkp)=members
!!$ else
!!$ ikpid(kpid)=0
!!$ endif
!!$ ! find coords of reduced onek
!!$ call coorskp(redkp(kpid),nkp)
!!$ enddo
!!$ enddo
!!$ enddo
allocate(done(div(1)*div(2)*div(3)))
done(:)=.false.
! debug information
if (tetradbglv.gt.0) then
write(*,*) 'libbzint(reduk): number of symmetries: ',nsymt
do i=1,nsymt
write(*,*) i
write(*,*) iio(1,:,i)
write(*,*) iio(2,:,i)
write(*,*) iio(3,:,i)
write(*,*)
end do
end if
do i1=0,div(1)-1
kk(1)=divsh*i1+shift(1)
onek(1)=i1
do i2=0,div(2)-1
kk(2)=divsh*i2+shift(2)
onek(2)=i2
do i3=0,div(3)-1
kk(3)=divsh*i3+shift(3)
onek(3)=i3
kpid=idkp(onek)
if (.not.done(kpid)) then
nirkp=nirkp+1
ikpid(kpid)=nirkp
do i=1,48
starr(i)=0
enddo
do i=1,nsymt
! from internal coordinates to lattice coordinates
kpr(:)=dble(kk(:))/dble(divsh*div(:))
! rotate k-point
kprt(:)=iio(:,1,i)*kpr(1)+iio(:,2,i)*kpr(2)+&
iio(:,3,i)*kpr(3)
! map to reciprocal unit cell
call mapto01(epslat,kprt)
! from lattice coordinates to internal coordinates
! for unshifted mesh
kprt(:)=(kprt(:)*div(:)*divsh-shift(:))/dble(divsh)
! location of rotated k-point on integer grid
nkp(:)=nint(kprt(:))
! determine fractional part of rotated k-point
call mapto01(epslat,kprt)
! if fractional part is present discard k-point
if (any(kprt.gt.epslat)) nkp(:)=-1
nkpid=idkp(nkp)
if (.not.all(nkp.eq.-1)) then
done(nkpid)=.true.
redkp(nkpid)=kpid
starr(i)=nkpid
end if
enddo
members=0
do i=1,nsymt
if(starr(i).ne.0) then
members=members+1
do j=i+1,nsymt
if(starr(i).eq.starr(j)) starr(j)=0
enddo
endif
enddo
weight(nirkp)=members
else
ikpid(kpid)=0
endif
! find coords of reduced onek
call coorskp(redkp(kpid),nkp)
enddo
enddo
enddo
deallocate(done)
else if (trim(tetraifc)=='exciting') then
! new code
allocate(done(div(1)*div(2)*div(3)))
done(:)=.false.
! debug information
if (tetradbglv.gt.0) then
write(*,*) 'libbzint(reduk): number of symmetries: ',nsymt
do i=1,nsymt
write(*,*) i
write(*,*) iio(1,:,i)
write(*,*) iio(2,:,i)
write(*,*) iio(3,:,i)
write(*,*)
end do
end if
do i3=0,div(3)-1
kk(3)=divsh*i3+shift(3)
onek(3)=i3
do i2=0,div(2)-1
kk(2)=divsh*i2+shift(2)
onek(2)=i2
do i1=0,div(1)-1
kk(1)=divsh*i1+shift(1)
onek(1)=i1
kpid=idkp(onek)
if (.not.done(kpid)) then
nirkp=nirkp+1
ikpid(kpid)=nirkp
do i=1,48
starr(i)=0
enddo
do i=1,nsymt
! from internal coordinates to lattice coordinates
kpr(:)=dble(kk(:))/dble(divsh*div(:))
! rotate k-point
kprt(:)=iio(:,1,i)*kpr(1)+iio(:,2,i)*kpr(2)+&
iio(:,3,i)*kpr(3)
! map to reciprocal unit cell
call mapto01(epslat,kprt)
! from lattice coordinates to internal coordinates
! for unshifted mesh
kprt(:)=(kprt(:)*div(:)*divsh-shift(:))/dble(divsh)
! location of rotated k-point on integer grid
nkp(:)=nint(kprt(:))
! determine fractional part of rotated k-point
call mapto01(epslat,kprt)
! if fractional part is present discard k-point
if (any(kprt.gt.epslat)) nkp(:)=-1
nkpid=idkp(nkp)
if (.not.all(nkp.eq.-1)) then
done(nkpid)=.true.
redkp(nkpid)=kpid
starr(i)=nkpid
end if
enddo
members=0
do i=1,nsymt
if(starr(i).ne.0) then
members=members+1
do j=i+1,nsymt
if(starr(i).eq.starr(j)) starr(j)=0
enddo
endif
enddo
weight(nirkp)=members
else
ikpid(kpid)=0
endif
! find coords of reduced onek
call coorskp(redkp(kpid),nkp)
enddo
enddo
enddo
deallocate(done)
! end new code
end if ! if (tetraifc)
!</sag>
! *** DEBUG ***
if (tetradbglv.gt.1) then
write(*,*) 'REDUK REPORTS: kpid,ikpid(kpid),redkp(kpid)'
write(*,*) 'div',div
do kpid=1,div(1)*div(2)*div(3)
write(*,*) kpid,ikpid(kpid),redkp(kpid)
end do
write(*,*)
end if
deallocate(iio)
! write(6,*)'------------------------------------------------------'
! write(6,*)' reduk: end'
! write(6,*)'------------------------------------------------------'
return
!<sag>
contains
subroutine mapto01(epsl,v)
implicit none
! arguments
real(8), intent(in) :: epsl
real(8), intent(inout) :: v(3)
! local variables
integer :: k
do k=1,3
v(k)=v(k)-dble(int(v(k)))
if (v(k).lt.0.d0) v(k)=v(k)+1.d0
if ((1.d0-v(k).lt.epsl).or.(v(k).lt.epsl)) v(k)=0.d0
end do
end subroutine mapto01
!</sag>
end subroutine reduk
!EOC