The Elk Code
genkpakq.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2017 T. Mueller, J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine genkpakq
7 use modmain
8 use modulr
9 implicit none
10 ! local variables
11 integer i1,i2,i3,j1,j2,j3
12 integer ikpa,ik0,ik
13 integer n1,iq,ifq,ir
14 real(8) v(3)
15 ! allocatable arrays
16 real(8), allocatable :: vkl0(:,:),vkc0(:,:),wkpt0(:)
17 ! find next largest FFT-compatible Q-point grid size
18 call nfftifc(npfftq,3,ngridq)
19 ! total number of Q-points
20 nqpt=ngridq(1)*ngridq(2)*ngridq(3)
21 ! number of complex FFT elements for real-complex transforms
22 n1=ngridq(1)/2+1
23 nfqrz=n1*ngridq(2)*ngridq(3)
24 ! integer grid intervals for the Q-points
25 intq(1,1:3)=ngridq(1:3)/2-ngridq(1:3)+1
26 intq(2,1:3)=ngridq(1:3)/2
27 ! κ-point grid should be half the Q-point grid
28 ngridkpa(1:3)=(ngridq(1:3)+1)/2
29 ! number of κ-points
31 ! integer grid intervals for the κ-points
32 intkpa(1,1:3)=ngridkpa(1:3)/2-ngridkpa(1:3)+1
33 intkpa(2,1:3)=ngridkpa(1:3)/2
34 ! allocate global Q-point arrays
35 if (allocated(ivq)) deallocate(ivq)
36 allocate(ivq(3,nqpt))
37 if (allocated(ivqiq)) deallocate(ivqiq)
38 allocate(ivqiq(intq(1,1):intq(2,1),intq(1,2):intq(2,2),intq(1,3):intq(2,3)))
39 if (allocated(iqfft)) deallocate(iqfft)
40 allocate(iqfft(nqpt))
41 if (allocated(ifqrz)) deallocate(ifqrz)
42 allocate(ifqrz(nqpt))
43 if (allocated(iqrzf)) deallocate(iqrzf)
44 allocate(iqrzf(nfqrz))
45 if (allocated(vql)) deallocate(vql)
46 allocate(vql(3,nqpt))
47 if (allocated(vqc)) deallocate(vqc)
48 allocate(vqc(3,nqpt))
49 ! store the κ-points as the first nkpa entries in the Q-point arrays
50 iq=0
51 do i1=intkpa(1,1),intkpa(2,1)
52  do i2=intkpa(1,2),intkpa(2,2)
53  do i3=intkpa(1,3),intkpa(2,3)
54  iq=iq+1
55  ivq(1,iq)=i1
56  ivq(2,iq)=i2
57  ivq(3,iq)=i3
58  end do
59  end do
60 end do
61 ! store the remaining Q-points
62 do i1=intq(1,1),intq(2,1)
63  do i2=intq(1,2),intq(2,2)
64  do i3=intq(1,3),intq(2,3)
65  if ((i1 < intkpa(1,1)).or.(i1 > intkpa(2,1)).or. &
66  (i2 < intkpa(1,2)).or.(i2 > intkpa(2,2)).or. &
67  (i3 < intkpa(1,3)).or.(i3 > intkpa(2,3))) then
68  iq=iq+1
69  ivq(1,iq)=i1
70  ivq(2,iq)=i2
71  ivq(3,iq)=i3
72  end if
73  end do
74  end do
75 end do
76 ! ensure the first point is the zero vector
77 do iq=1,nkpa
78  if ((ivq(1,iq) == 0).and.(ivq(2,iq) == 0).and.(ivq(3,iq) == 0)) then
79  ivq(1:3,iq)=ivq(1:3,1)
80  ivq(1:3,1)=0
81  exit
82  end if
83 end do
84 do iq=1,nqpt
85  i1=ivq(1,iq); i2=ivq(2,iq); i3=ivq(3,iq)
86 ! map from (i1,i2,i3) to Q-vector index
87  ivqiq(i1,i2,i3)=iq
88 ! Q-vector in Cartesian coordinates
89  vqc(1:3,iq)=dble(i1)*bvecu(1:3,1) &
90  +dble(i2)*bvecu(1:3,2) &
91  +dble(i3)*bvecu(1:3,3)
92 ! Q-vector in (unit cell) lattice coordinates
93  call r3mv(binv,vqc(:,iq),vql(:,iq))
94  where(abs(vql(1:3,iq)) < epslat) vql(1:3,iq)=0.d0
95 end do
96 ! set up Fourier transform index
97 do iq=1,nqpt
98  i1=ivq(1,iq); i2=ivq(2,iq); i3=ivq(3,iq)
99  if (i1 >= 0) then
100  j1=i1
101  else
102  j1=ngridq(1)+i1
103  end if
104  if (i2 >= 0) then
105  j2=i2
106  else
107  j2=ngridq(2)+i2
108  end if
109  if (i3 >= 0) then
110  j3=i3
111  else
112  j3=ngridq(3)+i3
113  end if
114  iqfft(iq)=j3*ngridq(2)*ngridq(1)+j2*ngridq(1)+j1+1
115 ! map from q-point index to real-complex FFT index and vice versa
116  if (i1 >= 0) then
117  ifq=j3*ngridq(2)*n1+j2*n1+j1+1
118  ifqrz(iq)=ifq
119  iqrzf(ifq)=iq
120  end if
121 end do
122 ! store the R-vectors in Cartesian coordinates spanning the ultracell
123 if (allocated(vrcu)) deallocate(vrcu)
124 allocate(vrcu(3,nqpt))
125 ir=0
126 do i3=0,ngridq(3)-1
127  v(3)=dble(i3)/dble(ngridq(3))
128  do i2=0,ngridq(2)-1
129  v(2)=dble(i2)/dble(ngridq(2))
130  do i1=0,ngridq(1)-1
131  v(1)=dble(i1)/dble(ngridq(1))
132  ir=ir+1
133  call r3mv(avecu,v,vrcu(:,ir))
134  end do
135  end do
136 end do
137 ! allocate the k-point weight array for band structure calculation
138 if (any(task == [720,725])) then
139  if (allocated(wkpt)) deallocate(wkpt)
140  allocate(wkpt(nkpt))
141  wkpt(1:nkpt)=1.d0
142 end if
143 ! store the existing k-point and weight arrays
144 allocate(vkl0(3,nkpt),vkc0(3,nkpt),wkpt0(nkpt))
145 vkl0(1:3,1:nkpt)=vkl(1:3,1:nkpt)
146 vkc0(1:3,1:nkpt)=vkc(1:3,1:nkpt)
147 wkpt0(1:nkpt)=wkpt(1:nkpt)
148 ! number of k+κ-points
149 nkpt0=nkpt
151 ! deallocate and reallocate k-point and weight arrays
152 deallocate(vkl,vkc,wkpt)
153 allocate(vkl(3,nkpt),vkc(3,nkpt),wkpt(nkpt))
154 ik=0
155 do ik0=1,nkpt0
156  do ikpa=1,nkpa
157  ik=ik+1
158  vkl(1:3,ik)=vkl0(1:3,ik0)+vql(1:3,ikpa)
159  vkc(1:3,ik)=vkc0(1:3,ik0)+vqc(1:3,ikpa)
160  wkpt(ik)=wkpt0(ik0)/dble(nkpa)
161  end do
162 end do
163 deallocate(vkl0,vkc0,wkpt0)
164 end subroutine
165 
integer task
Definition: modmain.f90:1299
integer nqpt
Definition: modmain.f90:525
integer nkpt
Definition: modmain.f90:461
integer, dimension(:,:), allocatable ivq
Definition: modmain.f90:529
integer nkpt0
Definition: modulr.f90:18
integer, dimension(3) ngridkpa
Definition: modulr.f90:20
integer, dimension(:), allocatable iqrzf
Definition: modmain.f90:543
real(8), dimension(3, 3) bvecu
Definition: modulr.f90:14
integer, dimension(:), allocatable iqfft
Definition: modmain.f90:537
real(8), dimension(:,:), allocatable vkc
Definition: modmain.f90:473
subroutine genkpakq
Definition: genkpakq.f90:7
real(8), dimension(:,:), allocatable vqc
Definition: modmain.f90:547
real(8), dimension(:), allocatable wkpt
Definition: modmain.f90:475
real(8), dimension(:,:), allocatable vql
Definition: modmain.f90:545
real(8), dimension(3, 3) avecu
Definition: modulr.f90:12
integer, dimension(:,:,:), allocatable ivqiq
Definition: modmain.f90:531
integer npfftq
Definition: modmain.f90:535
integer, dimension(3) ngridq
Definition: modmain.f90:515
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8), dimension(3, 3) binv
Definition: modmain.f90:18
integer nfqrz
Definition: modmain.f90:539
integer, dimension(2, 3) intq
Definition: modmain.f90:517
real(8) epslat
Definition: modmain.f90:24
real(8), dimension(:,:), allocatable vrcu
Definition: modulr.f90:26
integer, dimension(:), allocatable ifqrz
Definition: modmain.f90:541
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
integer, dimension(2, 3) intkpa
Definition: modulr.f90:22
Definition: modulr.f90:6
integer nkpa
Definition: modulr.f90:24
subroutine nfftifc(np, nd, n)
Definition: nfftifc.f90:10