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 ik0,ik,ikpa,nk,i
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 unless specified
28 do i=1,3
29  nk=(ngridq(i)+1)/2
30  if (ngridkpa(i) > 0) then
31  ngridkpa(i)=min(ngridkpa(i),nk)
32  else
33  ngridkpa(i)=nk
34  end if
35 end do
36 ! number of κ-points
38 ! integer grid intervals for the κ-points
39 intkpa(1,1:3)=ngridkpa(1:3)/2-ngridkpa(1:3)+1
40 intkpa(2,1:3)=ngridkpa(1:3)/2
41 ! allocate global Q-point arrays
42 if (allocated(ivq)) deallocate(ivq)
43 allocate(ivq(3,nqpt))
44 if (allocated(ivqiq)) deallocate(ivqiq)
45 allocate(ivqiq(intq(1,1):intq(2,1),intq(1,2):intq(2,2),intq(1,3):intq(2,3)))
46 if (allocated(iqfft)) deallocate(iqfft)
47 allocate(iqfft(nqpt))
48 if (allocated(ifqrz)) deallocate(ifqrz)
49 allocate(ifqrz(nqpt))
50 if (allocated(iqrzf)) deallocate(iqrzf)
51 allocate(iqrzf(nfqrz))
52 if (allocated(vql)) deallocate(vql)
53 allocate(vql(3,nqpt))
54 if (allocated(vqc)) deallocate(vqc)
55 allocate(vqc(3,nqpt))
56 ! store the κ-points as the first nkpa entries in the Q-point arrays
57 iq=0
58 do i1=intkpa(1,1),intkpa(2,1)
59  do i2=intkpa(1,2),intkpa(2,2)
60  do i3=intkpa(1,3),intkpa(2,3)
61  iq=iq+1
62  ivq(1,iq)=i1
63  ivq(2,iq)=i2
64  ivq(3,iq)=i3
65  end do
66  end do
67 end do
68 ! store the remaining Q-points
69 do i1=intq(1,1),intq(2,1)
70  do i2=intq(1,2),intq(2,2)
71  do i3=intq(1,3),intq(2,3)
72  if ((i1 < intkpa(1,1)).or.(i1 > intkpa(2,1)).or. &
73  (i2 < intkpa(1,2)).or.(i2 > intkpa(2,2)).or. &
74  (i3 < intkpa(1,3)).or.(i3 > intkpa(2,3))) then
75  iq=iq+1
76  ivq(1,iq)=i1
77  ivq(2,iq)=i2
78  ivq(3,iq)=i3
79  end if
80  end do
81  end do
82 end do
83 ! ensure the first point is the zero vector
84 do iq=1,nkpa
85  if ((ivq(1,iq) == 0).and.(ivq(2,iq) == 0).and.(ivq(3,iq) == 0)) then
86  ivq(1:3,iq)=ivq(1:3,1)
87  ivq(1:3,1)=0
88  exit
89  end if
90 end do
91 do iq=1,nqpt
92  i1=ivq(1,iq); i2=ivq(2,iq); i3=ivq(3,iq)
93 ! map from (i1,i2,i3) to Q-vector index
94  ivqiq(i1,i2,i3)=iq
95 ! Q-vector in Cartesian coordinates
96  vqc(1:3,iq)=dble(i1)*bvecu(1:3,1) &
97  +dble(i2)*bvecu(1:3,2) &
98  +dble(i3)*bvecu(1:3,3)
99 ! Q-vector in (unit cell) lattice coordinates
100  call r3mv(binv,vqc(:,iq),vql(:,iq))
101  where(abs(vql(1:3,iq)) < epslat) vql(1:3,iq)=0.d0
102 end do
103 ! set up Fourier transform index
104 do iq=1,nqpt
105  i1=ivq(1,iq); i2=ivq(2,iq); i3=ivq(3,iq)
106  if (i1 >= 0) then
107  j1=i1
108  else
109  j1=ngridq(1)+i1
110  end if
111  if (i2 >= 0) then
112  j2=i2
113  else
114  j2=ngridq(2)+i2
115  end if
116  if (i3 >= 0) then
117  j3=i3
118  else
119  j3=ngridq(3)+i3
120  end if
121  iqfft(iq)=j3*ngridq(2)*ngridq(1)+j2*ngridq(1)+j1+1
122 ! map from q-point index to real-complex FFT index and vice versa
123  if (i1 >= 0) then
124  ifq=j3*ngridq(2)*n1+j2*n1+j1+1
125  ifqrz(iq)=ifq
126  iqrzf(ifq)=iq
127  end if
128 end do
129 ! store the R-vectors in Cartesian coordinates spanning the ultracell
130 if (allocated(vrcu)) deallocate(vrcu)
131 allocate(vrcu(3,nqpt))
132 ir=0
133 do i3=0,ngridq(3)-1
134  v(3)=dble(i3)/dble(ngridq(3))
135  do i2=0,ngridq(2)-1
136  v(2)=dble(i2)/dble(ngridq(2))
137  do i1=0,ngridq(1)-1
138  v(1)=dble(i1)/dble(ngridq(1))
139  ir=ir+1
140  call r3mv(avecu,v,vrcu(:,ir))
141  end do
142  end do
143 end do
144 ! allocate the k-point weight array for band structure calculation
145 if (any(task == [720,725])) then
146  if (allocated(wkpt)) deallocate(wkpt)
147  allocate(wkpt(nkpt))
148  wkpt(1:nkpt)=1.d0
149 end if
150 ! store the existing k-point and weight arrays
151 allocate(vkl0(3,nkpt),vkc0(3,nkpt),wkpt0(nkpt))
152 vkl0(1:3,1:nkpt)=vkl(1:3,1:nkpt)
153 vkc0(1:3,1:nkpt)=vkc(1:3,1:nkpt)
154 wkpt0(1:nkpt)=wkpt(1:nkpt)
155 ! number of k+κ-points
156 nkpt0=nkpt
158 ! deallocate and reallocate k-point and weight arrays
159 deallocate(vkl,vkc,wkpt)
160 allocate(vkl(3,nkpt),vkc(3,nkpt),wkpt(nkpt))
161 ik=0
162 do ik0=1,nkpt0
163  do ikpa=1,nkpa
164  ik=ik+1
165  vkl(1:3,ik)=vkl0(1:3,ik0)+vql(1:3,ikpa)
166  vkc(1:3,ik)=vkc0(1:3,ik0)+vqc(1:3,ikpa)
167  wkpt(ik)=wkpt0(ik0)/dble(nkpa)
168  end do
169 end do
170 deallocate(vkl0,vkc0,wkpt0)
171 end subroutine
172 
integer task
Definition: modmain.f90:1297
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