TABLE OF CONTENTS


ABINIT/d2c_wtq [ Functions ]

[ Top ] [ Functions ]

NAME

  d2c_wtq

FUNCTION

  This routine calculates the integration weights on a coarse k-grid 
   using the integration weights from a denser k-gird. The weights of 
   the extra k points that being shared are evenly distributed.

COPYRIGHT

  Copyright (C) 2012-2018 ABINIT group (BXU)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

  elph_ds%k_fine%nkpt = number of fine q-points
  elph_ds%k_fine%wtq = integration weights of the fine q-grid
  elph_ds%k_phon%nkpt = number of coarse q-points

OUTPUT

  elph_ds%k_phon%wtq = integration weights of the coarse k-grid

PARENTS

      mka2f,mka2f_tr

CHILDREN

SOURCE

 32 #if defined HAVE_CONFIG_H
 33 #include "config.h"
 34 #endif
 35 
 36 #include "abi_common.h"
 37 
 38 subroutine d2c_wtq(elph_ds)
 39     
 40  use defs_basis
 41  use defs_elphon
 42  use m_profiling_abi
 43  use m_errors
 44 
 45 !This section has been created automatically by the script Abilint (TD).
 46 !Do not modify the following lines by hand.
 47 #undef ABI_FUNC
 48 #define ABI_FUNC 'd2c_wtq'
 49 !End of the abilint section
 50 
 51  implicit none
 52 
 53 !Arguments ------------------------------------
 54  type(elph_type),intent(inout) :: elph_ds
 55 
 56 !Local variables-------------------------------
 57  integer :: ii, jj, kk 
 58  integer :: ikpt, jkpt, kkpt
 59  integer :: iikpt, jjkpt, kkkpt
 60  integer :: ikpt_fine, ikpt_phon
 61  integer :: nkpt_fine1, nkpt_phon1
 62  integer :: nkpt_fine2, nkpt_phon2
 63  integer :: nkpt_fine3, nkpt_phon3
 64  integer :: nscale1, nscale2, nscale3
 65  
 66 ! *************************************************************************
 67  nkpt_phon1 = elph_ds%kptrlatt(1,1)
 68  nkpt_phon2 = elph_ds%kptrlatt(2,2)
 69  nkpt_phon3 = elph_ds%kptrlatt(3,3)
 70  nkpt_fine1 = elph_ds%kptrlatt_fine(1,1)
 71  nkpt_fine2 = elph_ds%kptrlatt_fine(2,2)
 72  nkpt_fine3 = elph_ds%kptrlatt_fine(3,3)
 73  nscale1 = dble(nkpt_fine1/nkpt_phon1)
 74  nscale2 = dble(nkpt_fine2/nkpt_phon2)
 75  nscale3 = dble(nkpt_fine3/nkpt_phon3)
 76  if (abs(INT(nscale1)-nscale1) > 0.01) then
 77    MSG_ERROR('The denser k-gird MUST be multiples of the phon k-grid')
 78  end if
 79  if (abs(INT(nscale2)-nscale2) > 0.01) then
 80    MSG_ERROR('The denser k-gird MUST be multiples of the phon k-grid')
 81  end if
 82  if (abs(INT(nscale3)-nscale3) > 0.01) then
 83    MSG_ERROR('The denser k-gird MUST be multiples of the phon k-grid')
 84  end if
 85  nscale1 = INT(nscale1)
 86  nscale2 = INT(nscale2)
 87  nscale3 = INT(nscale3)
 88  
 89 !bxu, get wtq of coarse grid from fine grid
 90  elph_ds%k_phon%wtq = zero
 91 
 92  do ikpt = 1, nkpt_phon1
 93    do jkpt = 1, nkpt_phon2
 94      do kkpt = 1, nkpt_phon3
 95        ikpt_phon = kkpt + (jkpt-1)*nkpt_phon3 + (ikpt-1)*nkpt_phon2*nkpt_phon3
 96 !      inside the paralellepipe
 97        do ii = -((nscale1+1)/2-1), ((nscale1+1)/2-1)
 98          do jj = -((nscale2+1)/2-1), ((nscale2+1)/2-1)
 99            do kk = -((nscale3+1)/2-1), ((nscale3+1)/2-1)
100              iikpt = 1 + (ikpt-1)*nscale1 + ii
101              jjkpt = 1 + (jkpt-1)*nscale2 + jj
102              kkkpt = 1 + (kkpt-1)*nscale3 + kk
103              if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
104              if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
105              if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
106              if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
107              if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
108              if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
109              ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
110              elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
111 &             elph_ds%k_fine%wtq(:,ikpt_fine,:)
112            end do
113          end do
114        end do
115 !      on the 6 faces
116        if (MOD(nscale3,2) == 0) then ! when nscale3 is an even number
117          do ii = -((nscale1+1)/2-1), ((nscale1+1)/2-1)
118            do jj = -((nscale2+1)/2-1), ((nscale2+1)/2-1)
119              iikpt = 1 + (ikpt-1)*nscale1 + ii
120              jjkpt = 1 + (jkpt-1)*nscale2 + jj
121              if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
122              if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
123              if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
124              if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
125 
126              kkkpt = 1 + (kkpt-1)*nscale3 + nscale3/2
127              if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
128              if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
129              ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
130              elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
131 &             0.5_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
132 
133              kkkpt = 1 + (kkpt-1)*nscale3 - nscale3/2
134              if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
135              if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
136              ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
137              elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
138 &             0.5_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
139            end do
140          end do
141        end if 
142        if (MOD(nscale2,2) == 0) then ! when nscale2 is an even number
143          do ii = -((nscale1+1)/2-1), ((nscale1+1)/2-1)
144            do kk = -((nscale3+1)/2-1), ((nscale3+1)/2-1)
145              iikpt = 1 + (ikpt-1)*nscale1 + ii
146              kkkpt = 1 + (kkpt-1)*nscale3 + kk
147              if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
148              if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
149              if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
150              if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
151 
152              jjkpt = 1 + (jkpt-1)*nscale2 + nscale2/2
153              if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
154              if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
155              ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
156              elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
157 &             0.5_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
158 
159              jjkpt = 1 + (jkpt-1)*nscale2 - nscale2/2
160              if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
161              if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
162              ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
163              elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
164 &             0.5_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
165            end do
166          end do
167        end if 
168        if (MOD(nscale1,2) == 0) then ! when nscale1 is an even number
169          do kk = -((nscale3+1)/2-1), ((nscale3+1)/2-1)
170            do jj = -((nscale2+1)/2-1), ((nscale2+1)/2-1)
171              kkkpt = 1 + (kkpt-1)*nscale3 + kk
172              jjkpt = 1 + (jkpt-1)*nscale2 + jj
173              if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
174              if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
175              if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
176              if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
177 
178              iikpt = 1 + (ikpt-1)*nscale1 + nscale1/2
179              if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
180              if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
181              ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
182              elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
183 &             0.5_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
184 
185              iikpt = 1 + (ikpt-1)*nscale1 - nscale1/2
186              if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
187              if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
188              ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
189              elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
190 &             0.5_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
191            end do
192          end do
193 !        on the 12 sides
194        end if 
195        if (MOD(nscale2,2) == 0 .and. MOD(nscale3,2) == 0) then
196          do ii = -((nscale1+1)/2-1), ((nscale1+1)/2-1)
197            iikpt = 1 + (ikpt-1)*nscale1 + ii
198            if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
199            if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
200 
201            jjkpt = 1 + (jkpt-1)*nscale2 + nscale2/2
202            kkkpt = 1 + (kkpt-1)*nscale3 + nscale3/2
203            if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
204            if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
205            if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
206            if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
207            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
208            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
209 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
210 
211            jjkpt = 1 + (jkpt-1)*nscale2 + nscale2/2
212            kkkpt = 1 + (kkpt-1)*nscale3 - nscale3/2
213            if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
214            if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
215            if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
216            if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
217            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
218            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
219 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
220 
221            jjkpt = 1 + (jkpt-1)*nscale2 - nscale2/2
222            kkkpt = 1 + (kkpt-1)*nscale3 + nscale3/2
223            if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
224            if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
225            if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
226            if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
227            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
228            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
229 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
230 
231            jjkpt = 1 + (jkpt-1)*nscale2 - nscale2/2
232            kkkpt = 1 + (kkpt-1)*nscale3 - nscale3/2
233            if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
234            if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
235            if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
236            if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
237            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
238            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
239 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
240          end do
241        end if 
242        if (MOD(nscale1,2) == 0 .and. MOD(nscale3,2) == 0) then
243          do jj = -((nscale2+1)/2-1), ((nscale2+1)/2-1)
244            jjkpt = 1 + (jkpt-1)*nscale2 + jj
245            if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
246            if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
247 
248            iikpt = 1 + (ikpt-1)*nscale1 + nscale1/2
249            kkkpt = 1 + (kkpt-1)*nscale3 + nscale3/2
250            if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
251            if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
252            if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
253            if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
254            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
255            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
256 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
257 
258            iikpt = 1 + (ikpt-1)*nscale1 + nscale1/2
259            kkkpt = 1 + (kkpt-1)*nscale3 - nscale3/2
260            if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
261            if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
262            if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
263            if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
264            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
265            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
266 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
267 
268            iikpt = 1 + (ikpt-1)*nscale1 - nscale1/2
269            kkkpt = 1 + (kkpt-1)*nscale3 + nscale3/2
270            if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
271            if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
272            if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
273            if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
274            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
275            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
276 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
277 
278            iikpt = 1 + (ikpt-1)*nscale1 - nscale1/2
279            kkkpt = 1 + (kkpt-1)*nscale3 - nscale3/2
280            if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
281            if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
282            if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
283            if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
284            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
285            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
286 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
287          end do
288        end if 
289        if (MOD(nscale2,2) == 0 .and. MOD(nscale1,2) == 0) then
290          do kk = -((nscale3+1)/2-1), ((nscale3+1)/2-1)
291            kkkpt = 1 + (kkpt-1)*nscale3 + kk
292            if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
293            if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
294 
295            jjkpt = 1 + (jkpt-1)*nscale2 + nscale2/2
296            iikpt = 1 + (ikpt-1)*nscale1 + nscale1/2
297            if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
298            if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
299            if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
300            if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
301            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
302            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
303 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
304 
305            jjkpt = 1 + (jkpt-1)*nscale2 + nscale2/2
306            iikpt = 1 + (ikpt-1)*nscale1 - nscale1/2
307            if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
308            if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
309            if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
310            if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
311            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
312            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
313 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
314 
315            jjkpt = 1 + (jkpt-1)*nscale2 - nscale2/2
316            iikpt = 1 + (ikpt-1)*nscale1 + nscale1/2
317            if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
318            if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
319            if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
320            if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
321            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
322            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
323 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
324 
325            jjkpt = 1 + (jkpt-1)*nscale2 - nscale2/2
326            iikpt = 1 + (ikpt-1)*nscale1 - nscale1/2
327            if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
328            if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
329            if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
330            if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
331            ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
332            elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
333 &           0.25_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
334          end do
335 !        on the 8 corners
336        end if 
337        if (MOD(nscale1,2) == 0 .and. MOD(nscale2,2) == 0 .and. MOD(nscale3,2) == 0) then
338          iikpt = 1 + (ikpt-1)*nscale1 + nscale1/2
339          jjkpt = 1 + (jkpt-1)*nscale2 + nscale2/2
340          kkkpt = 1 + (kkpt-1)*nscale3 + nscale3/2
341          if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
342          if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
343          if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
344          if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
345          if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
346          if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
347          ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
348          elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
349 &         0.125_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
350 
351          iikpt = 1 + (ikpt-1)*nscale1 + nscale1/2
352          jjkpt = 1 + (jkpt-1)*nscale2 + nscale2/2
353          kkkpt = 1 + (kkpt-1)*nscale3 - nscale3/2
354          if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
355          if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
356          if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
357          if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
358          if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
359          if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
360          ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
361          elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
362 &         0.125_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
363 
364          iikpt = 1 + (ikpt-1)*nscale1 + nscale1/2
365          jjkpt = 1 + (jkpt-1)*nscale2 - nscale2/2
366          kkkpt = 1 + (kkpt-1)*nscale3 + nscale3/2
367          if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
368          if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
369          if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
370          if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
371          if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
372          if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
373          ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
374          elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
375 &         0.125_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
376 
377          iikpt = 1 + (ikpt-1)*nscale1 + nscale1/2
378          jjkpt = 1 + (jkpt-1)*nscale2 - nscale2/2
379          kkkpt = 1 + (kkpt-1)*nscale3 - nscale3/2
380          if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
381          if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
382          if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
383          if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
384          if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
385          if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
386          ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
387          elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
388 &         0.125_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
389 
390          iikpt = 1 + (ikpt-1)*nscale1 - nscale1/2
391          jjkpt = 1 + (jkpt-1)*nscale2 + nscale2/2
392          kkkpt = 1 + (kkpt-1)*nscale3 + nscale3/2
393          if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
394          if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
395          if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
396          if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
397          if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
398          if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
399          ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
400          elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
401 &         0.125_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
402 
403          iikpt = 1 + (ikpt-1)*nscale1 - nscale1/2
404          jjkpt = 1 + (jkpt-1)*nscale2 + nscale2/2
405          kkkpt = 1 + (kkpt-1)*nscale3 - nscale3/2
406          if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
407          if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
408          if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
409          if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
410          if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
411          if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
412          ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
413          elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
414 &         0.125_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
415 
416          iikpt = 1 + (ikpt-1)*nscale1 - nscale1/2
417          jjkpt = 1 + (jkpt-1)*nscale2 - nscale2/2
418          kkkpt = 1 + (kkpt-1)*nscale3 + nscale3/2
419          if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
420          if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
421          if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
422          if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
423          if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
424          if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
425          ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
426          elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
427 &         0.125_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
428 
429          iikpt = 1 + (ikpt-1)*nscale1 - nscale1/2
430          jjkpt = 1 + (jkpt-1)*nscale2 - nscale2/2
431          kkkpt = 1 + (kkpt-1)*nscale3 - nscale3/2
432          if (iikpt .le. 0) iikpt = iikpt + nkpt_fine1
433          if (jjkpt .le. 0) jjkpt = jjkpt + nkpt_fine2
434          if (kkkpt .le. 0) kkkpt = kkkpt + nkpt_fine3
435          if (iikpt .gt. nkpt_fine1) iikpt = iikpt - nkpt_fine1
436          if (jjkpt .gt. nkpt_fine2) jjkpt = jjkpt - nkpt_fine2
437          if (kkkpt .gt. nkpt_fine3) kkkpt = kkkpt - nkpt_fine3
438          ikpt_fine = kkkpt + (jjkpt-1)*nkpt_fine3 + (iikpt-1)*nkpt_fine2*nkpt_fine3
439          elph_ds%k_phon%wtq(:,ikpt_phon,:)=elph_ds%k_phon%wtq(:,ikpt_phon,:)+&
440 &         0.125_dp*elph_ds%k_fine%wtq(:,ikpt_fine,:)
441        end if
442      end do
443    end do
444  end do
445 
446 !bxu, divide by nscale^3 to be consistent with the normalization of kpt_phon
447  elph_ds%k_phon%wtq = elph_ds%k_phon%wtq/nscale1/nscale2/nscale3
448 
449 end subroutine d2c_wtq