fft_interpolation.f90 7.13 KB
Newer Older
Marco Govoni's avatar
Marco Govoni committed
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
!
! Copyright (C) 2015-2016 M. Govoni 
! This file is distributed under the terms of the
! GNU General Public License. See the file `License'
! in the root directory of the present distribution,
! or http://www.gnu.org/copyleft/gpl.txt .
!
! This file is part of WEST.
!
! Contributors to this file: 
! He Ma, Marco Govoni
!
MODULE fourier_interpolation
 ! -----------------------------------------------------------------
 !
 IMPLICIT NONE
 !
 CONTAINS
 !
 SUBROUTINE get_G2R_mapping (n1, n2, n3, n, nl, ndim)
    ! 
    ! INPUT  : n1, n2, n3 = dimension of arbitrary R grid 
    !          n          = actual number of PW
    !          ndim       - dimensions of nl, 1: (n), 2:(n,-n)
    ! OUTPUT : nl         = computed mapping from G to R space (i,e. from [1,n] to [1, n1*n2*n3] )
    !
    USE kinds,         ONLY : DP
    USE gvect,         ONLY : g
    USE cell_base,     ONLY : at, bg, tpiba
    !
    IMPLICIT NONE
    ! 
    ! I/O
    !
    INTEGER, INTENT(IN)  :: n1, n2, n3, n, ndim
    INTEGER, INTENT(OUT) :: nl(n,ndim)        ! mapping from [1, npw] to [1, n1*n2*n3]  
    !
    ! Workspace
    !
    INTEGER  :: ig
    INTEGER  :: h, k, l, hmax, kmax, lmax, hidx, kidx, lidx
    !
    SELECT CASE(ndim) 
    CASE(1,2) 
      ! ok
    CASE DEFAULT 
      CALL errore("","Bad mapping given: only 1 and 2",1) 
    END SELECT
    !
    nl = 0
    !
    ! Maximum Miller indexes associated to the R grid 
    !
    hmax = (n1 - 1) / 2
    kmax = (n2 - 1) / 2
    lmax = (n3 - 1) / 2
    !
    DO ig = 1, n
       !
       ! Get the current Miller index
       !
       h = NINT(SUM(g(:,ig) * at(:,1)))
       k = NINT(SUM(g(:,ig) * at(:,2)))
       l = NINT(SUM(g(:,ig) * at(:,3)))
       !
       IF( ABS(h) <= hmax .AND. ABS(k) <= kmax .AND. ABS(l) <= lmax ) THEN 
          !
          ! derive position of G vector in (n1, n2, n3) grid
          !
          IF ( h >= 0 ) THEN
             hidx = h + 1
          ELSE
             hidx = h + n1 + 1
          ENDIF
          !
          IF ( k >= 0 ) THEN
             kidx = k + 1
          ELSE
             kidx = k + n2 + 1
          ENDIF
          !
          IF ( l >= 0 ) THEN
             lidx = l + 1
          ELSE
             lidx = l + n3 + 1
          ENDIF
          !
          nl(ig,1) = hidx + (kidx - 1) * n1 + (lidx - 1) * n1 * n2
          !
       ENDIF
       !
    ENDDO
    !
    IF (ndim == 2 ) THEN 
       DO ig = 1, n
          !
          ! Get the current Miller index
          !
          h = NINT(SUM(-g(:,ig) * at(:,1)))
          k = NINT(SUM(-g(:,ig) * at(:,2)))
          l = NINT(SUM(-g(:,ig) * at(:,3)))
          !
          IF( ABS(h) <= hmax .AND. ABS(k) <= kmax .AND. ABS(l) <= lmax ) THEN 
             !
             ! derive position of G vector in (n1, n2, n3) grid
             !
             IF ( h >= 0 ) THEN
                hidx = h + 1
             ELSE
                hidx = h + n1 + 1
             ENDIF
             !
             IF ( k >= 0 ) THEN
                kidx = k + 1
             ELSE
                kidx = k + n2 + 1
             ENDIF
             !
             IF ( l >= 0 ) THEN
                lidx = l + 1
             ELSE
                lidx = l + n3 + 1
             ENDIF
             !
             nl(ig,2) = hidx + (kidx - 1) * n1 + (lidx - 1) * n1 * n2
             !
          ENDIF
          !
       ENDDO
    ENDIF
    !
 END SUBROUTINE
 !
 !
 SUBROUTINE single_fwfft_fromArbitraryRGrid (fr, n1, n2, n3, n, nx, ndim, nl, fg, igk)
   !
   ! Note that the interpolation needs to be called by all processors within one band group 
   !
   ! FWFFT : R ---> G
   !
   ! INPUT  : n1, n2, n3 = dimension of arbitrary R grid 
   !          n     = actual number of PW
   !          nx    = leading dimendion for fg
   !          ndim  = 1, 2
   !          fr    = ONE COMPLEX array containing ONE function in R space (note that the array is not distributed, i.e. dimension = n1*n2*n3 )
   !          nl    = pre-computed mapping from G to R space (i,e. from [1,n] to [1, n1*n2*n3] )
   ! OUTPUT : fg    = ONE COMPLEX array containing ONE functions in G space (note that the array is distributed ) 
   !
   USE kinds,         ONLY : DP
   USE fft_scalar,    ONLY : cfft3d
   USE mp_bands,      ONLY : intra_bgrp_comm, me_bgrp
   USE mp,            ONLY : mp_bcast
   !
   ! I/O
   !
   INTEGER,     INTENT(IN)       :: n1, n2, n3, n, nx, ndim
   INTEGER,     INTENT(IN)       :: nl(n,ndim)
   COMPLEX(DP), INTENT(IN)       :: fr(n1*n2*n3)
   COMPLEX(DP), INTENT(OUT)      :: fg(nx)
   INTEGER, INTENT(IN), OPTIONAL :: igk(n)
   !
   ! Workspace
   !
   INTEGER :: ig, idx
   !
   ! FFT is serial...
   !
   IF ( me_bgrp == 0 ) THEN
      !
      CALL cfft3d( fr, n1, n2, n3, n1, n2, n3, 1, -1)
      !
   ENDIF
   !
   ! ... result is broadcast
   !
   CALL mp_bcast( fr, 0, intra_bgrp_comm )
   !
   IF( PRESENT(igk) ) THEN
      DO ig=1, n
         !
         idx = nl(igk(ig),1)
         IF (idx > 0) fg(ig) = fr(nl(igk(ig),1))
         !
      ENDDO
   ELSE
      DO ig=1, n
         !
         idx = nl(ig,1)
         IF (idx > 0) fg(ig) = fr(nl(ig,1))
         !
      ENDDO
   ENDIF
   !
   DO ig = (n+1), nx
      !
      fg(ig) = (0.0_DP,0.0_DP)
      !
   ENDDO
   !
 ENDSUBROUTINE
 !
 !
 SUBROUTINE single_invfft_toArbitraryRGrid (fr, n1, n2, n3, n, nx, ndim, nl, fg, igk)
   !
   ! Note that the interpolation needs to be called by all processors within one band group 
   !
   ! FWFFT : R ---> G
   !
   ! INPUT  : n1, n2, n3 = dimension of arbitrary R grid 
   !          n     = actual number of PW
   !          nx    = leading dimendion for fg
   !          ndim  = 1,2
   !          fr    = ONE COMPLEX array containing ONE function in R space (note that the array is not distributed, i.e. dimension = n1*n2*n3 )
   !          nl    = pre-computed mapping from G to R space (i,e. from [1,n] to [1, n1*n2*n3] )
   ! OUTPUT : fg    = ONE COMPLEX array containing ONE functions in G space (note that the array is distributed ) 
   !
   USE kinds,         ONLY : DP
   USE fft_scalar,    ONLY : cfft3d
   USE mp_bands,      ONLY : intra_bgrp_comm, me_bgrp
   USE mp,            ONLY : mp_bcast
   !
   ! I/O
   !
   INTEGER,     INTENT(IN) :: n1, n2, n3, n, nx, ndim
   INTEGER,     INTENT(IN) :: nl(n,ndim)
   COMPLEX(DP), INTENT(OUT) :: fr(n1*n2*n3)
   COMPLEX(DP), INTENT(IN):: fg(nx)
   INTEGER, INTENT(IN), OPTIONAL :: igk(n)
   !
   ! Workspace
   !
   INTEGER :: ig, idx
   !
   fr = 0._DP
   !
   IF( PRESENT(igk) ) THEN
      DO ig=1, n
         !
         idx = nl(igk(ig),1)
         IF (idx > 0) THEN 
            fr(nl(igk(ig),1)) = fg(ig)
         ENDIF
         !
      ENDDO
   ELSE 
      DO ig=1, n
         !
         idx = nl(ig,1)
         IF (idx > 0) THEN 
            fr(nl(ig,1)) = fg(ig)
         ENDIF
         !
      ENDDO
      IF ( ndim == 2 ) THEN
         DO ig=1, n
            !
            idx = nl(ig,2)
            IF (idx > 0) THEN 
               fr(nl(ig,2)) = DCONJG(fg(ig))
            ENDIF
            !
         ENDDO
      ENDIF
   ENDIF
   !
   ! ... result is gathered
   !
   CALL mp_sum( fr, intra_bgrp_comm )
   !
   ! FFT is serial...
   !
   IF ( me_bgrp == 0 ) THEN
      !
      CALL cfft3d( fr, n1, n2, n3, n1, n2, n3, 1, 1)
      !
   ENDIF
   !
 ENDSUBROUTINE
 !
END MODULE