Octopus
pfft.F90
Go to the documentation of this file.
1!! Copyright (C) 2011 J. Alberdi-Rodriguez, P. Garcia RisueƱo, M. Oliveira
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
23 use, intrinsic :: iso_c_binding
25 implicit none
26
27 private
28
29#ifdef HAVE_PFFT
30
31 ! make public some symbols from the pfft library
32 public :: &
33 pfft_cleanup, &
34 pfft_create_procmesh_2d, &
35 pfft_destroy_plan, &
36 pfft_execute, &
37 pfft_init, &
38 pfft_local_size_dft_r2c_3d, &
39 pfft_local_size_dft_3d, &
40 pfft_plan_dft_c2r_3d, &
41 pfft_plan_dft_r2c_3d, &
42 pfft_plan_dft_3d, &
43 pfft_int, &
44 pfft_ptrdiff_t, &
45 pfft_float, &
46 pfft_double, &
47 pfft_ldouble, &
48 pfft_unsigned, &
49 pfft_transposed_in, &
50 pfft_transposed_out, &
51 pfft_estimate, &
52 pfft_tune, &
53 pfft_destroy_input, &
54 pfft_measure, &
55 pfft_forward, &
56 pfft_backward
57
61#ifndef HAVE_FFTW3_MPI
62 integer(C_INTPTR_T), parameter :: FFTW_MPI_DEFAULT_BLOCK = 0
63#endif
64 include "pfft.f03"
65#endif
66
67end module pfft_params_oct_m
68
71module pfft_oct_m
72 use debug_oct_m
73 use global_oct_m
74 use, intrinsic :: iso_c_binding
75 use math_oct_m
78
79 implicit none
80
81 private
82
83 public :: &
89
90contains
91
92 ! ---------------------------------------------------------
94 subroutine pfft_decompose(n_proc, dim1, dim2)
95 integer, intent(in) :: n_proc
96 integer, intent(out) :: dim1
97 integer, intent(out) :: dim2
98
99 integer :: np, i
100
101 push_sub(pfft_decompose)
102
103 assert(n_proc > 0)
104
105 dim1 = 1
106 dim2 = 1
107 np = n_proc
108 i = n_proc-1
109
110 if (is_prime(n_proc)) then
111 dim1 = n_proc
112 else
113 do
114 if (i <= 1) exit
115 if (mod(np,i) /= 0.or.(.not. is_prime(i))) then
116 i=i-1
117 cycle
118 end if
119 np=np/i
120 if (dim1 <= dim2) then
121 dim1 = dim1*i
122 else
123 dim2 = dim2*i
124 end if
125 end do
126 end if
127
128 assert(dim1*dim2 == n_proc)
129
130 pop_sub(pfft_decompose)
131 end subroutine pfft_decompose
132
133 ! ---------------------------------------------------------
135 subroutine pfft_prepare_plan_r2c(plan, n, in, out, sign, flags, mpi_comm)
136 type(C_PTR), intent(out) :: plan
137 integer, intent(in) :: n(:)
138 real(real64), pointer, intent(inout) :: in(:,:,:)
139 complex(real64), pointer, intent(inout) :: out(:,:,:)
140 integer, intent(in) :: sign
141 integer, intent(in) :: flags
143 integer, intent(in) :: mpi_comm
144
145 integer(C_INTPTR_T) :: tmp_n(3)
146
147 push_sub(pfft_prepare_plan_r2c)
148 call profiling_in("PFFT_PLAN_R2C")
149
150#ifdef HAVE_PFFT
151 assert(sign == pfft_forward)
152#endif
153
154 tmp_n(1:3) = n(3:1:-1)
156#ifdef HAVE_PFFT
157 plan = pfft_plan_dft_r2c_3d(tmp_n,in,out,mpi_comm,sign,&
158 pfft_transposed_out + pfft_measure + pfft_destroy_input + pfft_tune)
159#else
160 plan = c_null_ptr
161#endif
162
163 call profiling_out("PFFT_PLAN_R2C")
165 end subroutine pfft_prepare_plan_r2c
166
167 ! ---------------------------------------------------------
169 subroutine pfft_prepare_plan_c2r(plan, n, in, out, sign, flags, mpi_comm)
170 type(C_PTR), intent(out) :: plan
171 integer, intent(in) :: n(:)
172 complex(real64), pointer, intent(inout) :: in(:,:,:)
173 real(real64), pointer, intent(inout) :: out(:,:,:)
174 integer, intent(in) :: sign
175 integer, intent(in) :: flags
177 integer, intent(in) :: mpi_comm
178
179 integer(C_INTPTR_T) :: tmp_n(3)
180 push_sub(pfft_prepare_plan_c2r)
181 call profiling_in("PFFT_PLAN_C2R")
182
183#ifdef HAVE_PFFT
184 assert(sign == pfft_backward)
185#endif
186
187 tmp_n(1:3) = n(3:1:-1)
188
189#ifdef HAVE_PFFT
190 plan = pfft_plan_dft_c2r_3d(tmp_n,in,out,mpi_comm, sign, &
191 pfft_transposed_in + pfft_measure + pfft_destroy_input + pfft_tune)
192 !call PDFFT(plan_dft_c2r_3d) (plan, tmp_n(1), in(1,1,1), out(1,1,1), mpi_comm, sign, &
193 ! PFFT_TRANSPOSED_IN + PFFT_MEASURE + PFFT_DESTROY_INPUT + PFFT_TUNE)
194#else
195 plan = c_null_ptr
196#endif
197
198 call profiling_out("PFFT_PLAN_C2R")
199 pop_sub(pfft_prepare_plan_c2r)
200 end subroutine pfft_prepare_plan_c2r
201
202 ! ---------------------------------------------------------
204 subroutine pfft_prepare_plan_c2c(plan, n, in, out, sign, flags, mpi_comm)
205 type(c_ptr), intent(out) :: plan
206 integer, intent(in) :: n(:)
207 complex(real64), pointer, intent(inout) :: in(:,:,:)
208 complex(real64), pointer, intent(inout) :: out(:,:,:)
209 integer, intent(in) :: sign
210 ! 1 Has to be FFTW_FORWARD or FFTW_BACKWARD
211 integer, intent(in) :: flags
213 integer, intent(in) :: mpi_comm
214
215#ifdef HAVE_PFFT
216 integer(C_INT) :: pfft_flags
217#endif
218 integer(C_INTPTR_T) :: tmp_n(3)
219
220 push_sub(pfft_prepare_plan_c2c)
221
222 tmp_n(1:3) = n(3:1:-1)
223
224#ifdef HAVE_PFFT
225 if (sign == pfft_forward) then
226 pfft_flags = pfft_transposed_out
227 else
228 pfft_flags = pfft_transposed_in
229 end if
230
231 plan = pfft_plan_dft_3d(tmp_n,in,out,mpi_comm, sign, pfft_flags)
232 !call PDFFT(plan_dft_3d) (plan, tmp_n(1), in(1,1,1), out(1,1,1), mpi_comm, sign, pfft_flags)
233#else
234 plan = c_null_ptr
235#endif
236
237 pop_sub(pfft_prepare_plan_c2c)
238 end subroutine pfft_prepare_plan_c2c
239
240 ! ---------------------------------------------------------
241 subroutine pfft_get_dims(rs_n_global, mpi_comm, is_real, alloc_size, fs_n_global, rs_n, fs_n, rs_istart, fs_istart)
242 integer, intent(in) :: rs_n_global(1:3)
243 integer, intent(in) :: mpi_comm
244 logical, intent(in) :: is_real
245 integer(int64), intent(out) :: alloc_size
246 integer, intent(out) :: fs_n_global(1:3)
247 integer, intent(out) :: rs_n(1:3)
248 integer, intent(out) :: fs_n(1:3)
249 integer, intent(out) :: rs_istart(1:3)
250 integer, intent(out) :: fs_istart(1:3)
251
252 integer(C_INTPTR_T) :: tmp_alloc_size, tmp_n(3)
253 integer(C_INTPTR_T) :: tmp_rs_n(3), tmp_rs_istart(3)
254 integer(C_INTPTR_T) :: tmp_fs_n(3), tmp_fs_istart(3)
255
256 push_sub(pfft_get_dims)
257
258 fs_n_global(1:3) = rs_n_global(1:3)
259 ! if calling the ISO_C_BINDING routines from PFFT, one must take into account that
260 ! the array order is the C array order
261 tmp_n(1:3) = rs_n_global(3:1:-1)
263 tmp_alloc_size = 0
264#ifdef HAVE_PFFT
265 if (is_real) then
266 tmp_alloc_size = pfft_local_size_dft_r2c_3d(tmp_n, mpi_comm, pfft_transposed_out, &
267 tmp_rs_n, tmp_rs_istart, tmp_fs_n, tmp_fs_istart)
268 !call PDFFT(local_size_dft_r2c_3d)(tmp_alloc_size, tmp_n(1), mpi_comm, PFFT_TRANSPOSED_OUT, &
269 ! tmp_rs_n(1), tmp_rs_istart(1), tmp_fs_n(1), tmp_fs_istart(1))
270 fs_n_global(1) = rs_n_global(1)/2 + 1
271 else
272 tmp_alloc_size = pfft_local_size_dft_3d(tmp_n,mpi_comm, pfft_transposed_out, &
273 tmp_rs_n,tmp_rs_istart,tmp_fs_n,tmp_fs_istart)
274 !call PDFFT(local_size_dft_3d)(tmp_alloc_size, tmp_n(1), mpi_comm, PFFT_TRANSPOSED_OUT, &
275 ! tmp_rs_n(1), tmp_rs_istart(1), tmp_fs_n(1), tmp_fs_istart(1))
276 end if
277#else
278 tmp_rs_n = 0
279 tmp_fs_n = 0
280 tmp_rs_istart = 0
281 tmp_fs_istart = 0
282#endif
283 alloc_size = int(tmp_alloc_size, int64)
284 rs_n(1:3) = int(tmp_rs_n(3:1:-1))
285 fs_n(1:3) = int(tmp_fs_n(3:1:-1))
286 rs_istart(1:3) = int(tmp_rs_istart(3:1:-1)+1)
287 fs_istart(1:3) = int(tmp_fs_istart(3:1:-1)+1)
288
289 pop_sub(pfft_get_dims)
290 end subroutine pfft_get_dims
291
292end module pfft_oct_m
293
294!! Local Variables:
295!! mode: f90
296!! coding: utf-8
297!! End:
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
logical function, public is_prime(n)
Definition: math.F90:992
The low level module to work with the PFFT library. http:
Definition: pfft.F90:164
subroutine, public pfft_prepare_plan_c2c(plan, n, in, out, sign, flags, mpi_comm)
Octopus subroutine to prepare a PFFT plan real to complex.
Definition: pfft.F90:298
subroutine, public pfft_prepare_plan_c2r(plan, n, in, out, sign, flags, mpi_comm)
Octopus subroutine to prepare a PFFT plan real to complex.
Definition: pfft.F90:263
subroutine, public pfft_decompose(n_proc, dim1, dim2)
Decompose all available processors in 2D processor grid, most equally possible.
Definition: pfft.F90:188
subroutine, public pfft_prepare_plan_r2c(plan, n, in, out, sign, flags, mpi_comm)
Octopus subroutine to prepare a PFFT plan real to complex.
Definition: pfft.F90:229
subroutine, public pfft_get_dims(rs_n_global, mpi_comm, is_real, alloc_size, fs_n_global, rs_n, fs_n, rs_istart, fs_istart)
Definition: pfft.F90:335
The includes for the PFFT.
Definition: pfft.F90:115
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552