FV3 Bundle
fft.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 
20 !> \file
21 !! \brief Contains the \ref fft_mod module
22 
23 
24 module fft_mod
25 
26 !###############################################################################
27 !> \defgroup fft_mod fft_mod
28 !!
29 !! \author Bruce Wyman <Bruce.Wyman@noaa.gov>
30 !!
31 !! \brief This module performs simultaneous fast Fourier transforms (FFTs)
32 !! between real grid space and complex Fourier space.
33 !!
34 !! This routine computes multiple 1-dimensional FFTs and inverse FFTs.
35 !! There are 2d and 3d versions between type real grid point space
36 !! and type complex Fourier space. There are single (32-bit) and
37 !! full (64-bit) versions.
38 !!
39 !! On Cray and SGI systems, vendor-specific scientific library
40 !! routines are used, otherwise a user may choose a NAG library version
41 !! or stand-alone version using Temperton's FFT.
42 !!
43 !! fft uses the SCILIB on SGICRAY, otherwise the NAG library or
44 !! a standalone version of Temperton's fft is used:
45 !! - SCFFTM and CSFFTM are used on Crays
46 !! - DZFFTM and ZDFFTM are used on SGIs
47 !!
48 !! The following NAG routines are used: c06fpf, c06gqf, c06fqf.
49 !! These routine names may be slightly different on different
50 !! platforms.
51 !!
52 !! <b> Modules Included: </b>
53 !!
54 !! <table>
55 !! <tr>
56 !! <th> Module Name </th>
57 !! <th> Conditions </th>
58 !! <th> Functions Included </th>
59 !! </tr>
60 !! <tr>
61 !! <td> platform_mod </td>
62 !! <td> </td>
63 !! <td> R8_KIND, R4_KIND </td>
64 !! </tr>
65 !! <tr>
66 !! <td> fms_mod </td>
67 !! <td> </td>
68 !! <td> write_version_number, error_mesg, FATAL <\td>
69 !! </tr>
70 !! <tr>
71 !! <td> fft99_mod </td>
72 !! <td> if not on SGICRAY </td>
73 !! <td> fft991, set99 </td>
74 !! </tr>
75 !! </table>
76 !!
77 
78 !-----------------------------------------------------------------------
79 !these are used to determine hardware/OS/compiler
80 
81 #ifdef __sgi
82 # ifdef _COMPILER_VERSION
83 !the MIPSPro compiler defines _COMPILER_VERSION
84 # define sgi_mipspro
85 # else
86 # define sgi_generic
87 # endif
88 #endif
89 
90 !fft uses the SCILIB on SGICRAY, and the NAG library otherwise
91 #if defined(_CRAY) || defined(sgi_mipspro)
92 # define SGICRAY
93 #endif
94 
95 use platform_mod, only: r8_kind, r4_kind
96 use fms_mod, only: write_version_number, &
97  error_mesg, fatal
98 #ifndef SGICRAY
99 #ifndef NAGFFT
100 use fft99_mod, only: fft991, set99
101 #endif
102 #endif
103 
104 implicit none
105 private
106 
107 !----------------- interfaces --------------------
108 
110 
111 
112 
113 !###############################################################################
114 !> \defgroup fft_grid_to_fourier fft_grid_to_fourier
115 !! \ingroup fft_mod
116 !!
117 !! \brief Given multiple sequences of real data values, this routine
118 !! computes the complex Fourier transform for all sequences.
119 !!
120 !! <b> Implemented by: </b>
121 !! - \ref fft_grid_to_fourier_float_2d
122 !! - \ref fft_grid_to_fourier_double_2d
123 !! - \ref fft_grid_to_fourier_float_3d
124 !! - \ref fft_grid_to_fourier_double_3d
125 !!
126 !! \code{.f90}
127 !! fourier = fft_grid_to_fourier ( grid )
128 !! \endcode
129 !!
130 !! \param [in] <grid> Multiple sequence of real data values. The first dimension
131 !! must be n+1 (where n is the size of a single sequence).
132 !!
133 !! \param [out] <fourier> Multiple sequences of transformed data in complex
134 !! Fourier space. The first dimension must equal n/2+1 (where n is the
135 !! size of a single sequence). The remaining dimensions must be the
136 !! same size as the input argument "grid".
137 !!
138 !! The complex Fourier components are passed in the following format (where n
139 !! is length of each real transform).
140 !! <table>
141 !! <tr>
142 !! <td> fourier (1) </td> <td> = cmplx ( a(0), b(0) ) </td>
143 !! </tr>
144 !! <tr>
145 !! <td> fourier (2) </td> <td> = cmplx ( a(1), b(1) ) </td>
146 !! </tr>
147 !! <tr>
148 !! <td> : </td> <td> : </td>
149 !! </tr>
150 !! <tr>
151 !! <td> fourier (n/2+1) </td> <td> = cmplx ( a(n/2), b(n/2) ) </td>
152 !! </tr>
153 !! </table>
154 !!
155 !! \throw FATAL `"fft_init must be called"` \n
156 !! The initialization routine fft_init must be called before routines
157 !! fft_grid_to_fourier.
158 !!
159 !! \throw FATAL `"size of first dimension of input data is wrong"` \n
160 !! The real grid point field must have a first dimension equal to n+1
161 !! (where n is the size of each real transform). This message occurs
162 !! when using the SGI/Cray fft.
163 !!
164 !! \throw FATAL `"length of input data too small"` \n
165 !! The real grid point field must have a first dimension equal to n
166 !! (where n is the size of each real transform). This message occurs
167 !! when using the NAG or Temperton fft.
168 !!
169 !! \throw FATAL `"float kind not supported for nag fft"` \n
170 !! 32-bit real data is not supported when using the NAG fft. You
171 !! may try modifying this part of the code by uncommenting the
172 !! calls to the NAG library or less consider using the Temperton fft.
173 !!
177 end interface
178 
179 
180 
181 !###############################################################################
182 !> \defgroup fft_fourier_to_grid fft_fourier_to_grid
183 !! \ingroup fft_mod
184 !!
185 !! \brief Given multiple sequences of Fourier space transforms,
186 !! this routine computes the inverse transform and returns
187 !! the real data values for all sequences.
188 !!
189 !! \code{.f90}
190 !! grid = fft_fourier_to_grid ( fourier )
191 !! \endcode
192 !!
193 !! \param [in] <fourier> Multiple sequence complex Fourier space transforms.
194 !! The first dimension must equal n/2+1 (where n is the
195 !! size of a single real data sequence).
196 !!
197 !! \param [out] <grid> Multiple sequence of real data values. The first
198 !! dimension must be n+1 (where n is the size of a single sequence).
199 !! The remaining dimensions must be the same size as the input
200 !! argument "fourier".
201 !!
202 !! \throw FATAL `"fft_init must be called"` \n
203 !! The initialization routine fft_init must be called before routines
204 !! fft_fourier_to_grid.
205 !!
206 !! \throw FATAL `"size of first dimension of input data is wrong"` \n
207 !! The complex Fourier field must have a first dimension equal to
208 !! n/2+1 (where n is the size of each real transform). This message
209 !! occurs when using the SGI/Cray fft.
210 !!
211 !! \throw FATAL `"length of input data too small"` \n
212 !! The complex Fourier field must have a first dimension greater
213 !! than or equal to n/2+1 (where n is the size of each real
214 !! transform). This message occurs when using the NAG or Temperton fft.
215 !!
216 !! \throw FATAL `"float kind not supported for nag fft"` \n
217 !! 32-bit real data is not supported when using the NAG fft. You may try
218 !! modifying this part of the code by uncommenting the calls to the NAG
219 !! library or less consider using the Temperton fft.
220 !!
224 end interface
225 
226 
227 
228 !---------------------- private data -----------------------------------
229 
230 ! tables for trigonometric constants and factors
231 ! (not all will be used)
232 real(R8_KIND), allocatable, dimension(:) :: table8
233 real(R4_KIND), allocatable, dimension(:) :: table4
234 real , allocatable, dimension(:) :: table99
235 integer , allocatable, dimension(:) :: ifax
236 
237 logical :: do_log =.true.
238 integer :: leng, leng1, leng2, lenc ! related to transform size
239 
240 logical :: module_is_initialized=.false.
241 
242 ! Include variable "version" to be written to log file.
243 #include<file_version.h>
244 
245 contains
246 
247 
248 
249 !###############################################################################
250 !> \fn fft_grid_to_fourier_float_2d
251 !! \ingroup fft_grid_to_fourier
252 !! \implements fft_grid_to_fourier
253 !!
254 !! \brief A function that returns multiple transforms in complex
255 !! fourier space, the first dimension must equal n/2+1 (where n is the
256 !! size of each real transform). The remaining dimensions must be the
257 !! same size as the input argument "grid".
258 !!
259 !! \code{.f90}
260 !! real (R4_KIND), intent(in), dimension(:,:) :: grid
261 !! complex(R4_KIND), dimension(lenc,size(grid,2)) :: fourier
262 !! \endcode
263 !!
264 !! \param [in] <grid> Multiple transforms in grid point space, the first
265 !! dimension must be n+1 (where n is the size of each real transform).
266 !!
267 !! \param [out] <fourier>
268 !!
269 !! \throw FATAL `"fft_init must be called"` \n
270 !! The initialization routine fft_init must be called before routines
271 !! fft_grid_to_fourier.
272 !!
273 !! \throw FATAL `"size of first dimension of input data is wrong"` \n
274 !! The real grid point field must have a first dimension equal to n+1
275 !! (where n is the size of each real transform). This message occurs
276 !! when using the SGI/Cray fft.
277 !!
278 !! \throw FATAL `"length of input data too small"` \n
279 !! The real grid point field must have a first dimension equal to n
280 !! (where n is the size of each real transform). This message occurs
281 !! when using the NAG or Temperton fft.
282 !!
283 !! \throw FATAL `"float kind not supported for nag fft"` \n
284 !! 32-bit real data is not supported when using the NAG fft. You
285 !! may try modifying this part of the code by uncommenting the
286 !! calls to the NAG library or less consider using the Temperton fft.
287 !!
288 function fft_grid_to_fourier_float_2d (grid) result (fourier)
290 !-----------------------------------------------------------------------
291 
292  real (R4_KIND), intent(in), dimension(:,:) :: grid
293  complex(R4_KIND), dimension(lenc,size(grid,2)) :: fourier
294 
295 #ifdef SGICRAY
296 # ifdef _CRAY
297 ! local storage for cray fft
298  real(R4_KIND), dimension((2*leng+4)*size(grid,2)) :: work
299 # else
300 ! local storage for sgi fft
301  real(R4_KIND), dimension(leng2) :: work
302 # endif
303 #else
304 # ifdef NAGFFT
305 ! local storage for nag fft
306  real(R4_KIND), dimension(size(grid,2),leng) :: data, work
307 # else
308 ! local storage for temperton fft
309  real, dimension(leng2,size(grid,2)) :: data
310  real, dimension(leng1,size(grid,2)) :: work
311 # endif
312 #endif
313 
314 #if defined(SGICRAY) || defined(NAGFFT)
315  real(R4_KIND) :: scale
316 #endif
317  integer :: j, k, num, len_grid
318 
319 !-----------------------------------------------------------------------
320 
321  if (.not.module_is_initialized) &
322  call error_handler ('fft_grid_to_fourier', &
323  'fft_init must be called.')
324 
325 !-----------------------------------------------------------------------
326 
327  len_grid = size(grid,1)
328 #ifdef SGICRAY
329  if (len_grid /= leng1) call error_handler ('fft_grid_to_fourier', &
330  'size of first dimension of input data is wrong')
331 #else
332  if (len_grid < leng) call error_handler ('fft_grid_to_fourier', &
333  'length of input data too small.')
334 #endif
335 !-----------------------------------------------------------------------
336 !----------------transform to fourier coefficients (+1)-----------------
337 
338  num = size(grid,2) ! number of transforms
339 
340 #ifdef SGICRAY
341 ! Cray/SGI fft
342  scale = 1./real(leng)
343 # ifdef _CRAY
344  call scfftm (-1,leng,num,scale, grid,leng1, fourier,lenc, &
345  table4, work, 0)
346 # else
347  call scfftm (-1,leng,num,scale, grid,leng1, fourier,lenc, &
348  table4, work, 0)
349 # endif
350 #else
351 # ifdef NAGFFT
352 ! NAG fft
353 ! will not allow float kind for NAG
354  call error_handler ('fft_grid_to_fourier', &
355  'float kind not supported for nag fft')
356  do j=1,size(grid,2)
357  data(j,1:leng) = grid(1:leng,j)
358  enddo
359  scale = 1./sqrt(float(leng))
360  data = data * scale
361  fourier(1,:) = cmplx( data(:,1), 0. )
362  do k=2,lenc-1
363  fourier(k,:) = cmplx( data(:,k), data(:,leng-k+2) )
364  enddo
365  fourier(lenc,:) = cmplx( data(:,lenc), 0. )
366 # else
367 ! Temperton fft
368  do j=1,num
369  data(1:leng,j) = grid(1:leng,j)
370  enddo
371  call fft991 (data,work,table99,ifax,1,leng2,leng,num,-1)
372  do j=1,size(grid,2)
373  do k=1,lenc
374  fourier(k,j) = cmplx( data(2*k-1,j), data(2*k,j) )
375  enddo
376  enddo
377 # endif
378 #endif
379 !-----------------------------------------------------------------------
380 
381  end function fft_grid_to_fourier_float_2d
382 
383 
384 
385 !###############################################################################
386 !> \fn fft_fourier_to_grid_float_2d
387 !! \implements fft_fourier_to_grid
388 !!
389 !! \brief A function that multiple transforms in grid point space, the first
390 !! dimension must be n+1 (where n is the size of each real transform).
391 !! The remaining dimensions must be the same size as the input
392 !! argument "fourier".
393 !!
394 !! \code{.f90}
395 !! complex(R4_KIND), intent(in), dimension(:,:) :: fourier
396 !! real (R4_KIND), dimension(leng1,size(fourier,2)) :: grid
397 !! \endcode
398 !!
399 !! \param [in] <fourier> Multiple transforms in complex fourier space, the
400 !! first dimension must equal n/2+1 (where n is the size of each
401 !! real transform).
402 !!
403 !! \param [out] <grid>
404 !!
405 !! \throw FATAL `"fft_init must be called"` \n
406 !! The initialization routine fft_init must be called before routines
407 !! fft_fourier_to_grid.
408 !!
409 !! \throw FATAL `"size of first dimension of input data is wrong"` \n
410 !! The complex Fourier field must have a first dimension equal to
411 !! n/2+1 (where n is the size of each real transform). This message
412 !! occurs when using the SGI/Cray fft.
413 !!
414 !! \throw FATAL `"length of input data too small"` \n
415 !! The complex Fourier field must have a first dimension greater
416 !! than or equal to n/2+1 (where n is the size of each real
417 !! transform). This message occurs when using the NAG or Temperton fft.
418 !!
419 !! \throw FATAL `"float kind not supported for nag fft"` \n
420 !! 32-bit real data is not supported when using the NAG fft. You may try
421 !! modifying this part of the code by uncommenting the calls to the NAG
422 !! library or less consider using the Temperton fft.
423 !!
424  function fft_fourier_to_grid_float_2d (fourier) result (grid)
426 !-----------------------------------------------------------------------
427 
428  complex(R4_KIND), intent(in), dimension(:,:) :: fourier
429  real (R4_KIND), dimension(leng1,size(fourier,2)) :: grid
430 
431 #ifdef SGICRAY
432 # ifdef _CRAY
433 ! local storage for cray fft
434  real(R4_KIND), dimension((2*leng+4)*size(fourier,2)) :: work
435 # else
436 ! local storage for sgi fft
437  real(R4_KIND), dimension(leng2) :: work
438 # endif
439 #else
440 # ifdef NAGFFT
441 ! local storage for nag fft
442  real(R4_KIND), dimension(size(fourier,2),leng) :: data, work
443 # else
444 ! local storage for temperton fft
445  real, dimension(leng2,size(fourier,2)) :: data
446  real, dimension(leng1,size(fourier,2)) :: work
447 # endif
448 #endif
449 
450 #if defined(SGICRAY) || defined(NAGFFT)
451  real(R4_KIND) :: scale
452 #endif
453  integer :: j, k, num, len_fourier
454 
455 !-----------------------------------------------------------------------
456 
457  if (.not.module_is_initialized) &
458  call error_handler ('fft_grid_to_fourier', &
459  'fft_init must be called.')
460 
461 !-----------------------------------------------------------------------
462 
463  len_fourier = size(fourier,1)
464  num = size(fourier,2) ! number of transforms
465 
466 #ifdef SGICRAY
467  if (len_fourier /= lenc) call error_handler ('fft_fourier_to_grid', &
468  'size of first dimension of input data is wrong')
469 #else
470  if (len_fourier < lenc) call error_handler ('fft_fourier_to_grid', &
471  'length of input data too small.')
472 #endif
473 !-----------------------------------------------------------------------
474 !----------------inverse transform to real space (-1)-------------------
475 
476 #ifdef SGICRAY
477 ! Cray/SGI fft
478  scale = 1.0
479 # ifdef _CRAY
480  call csfftm (+1,leng,num,scale, fourier,len_fourier, &
481  grid,leng1, table4, work, 0)
482 # else
483  call csfftm (+1,leng,num,scale, fourier,len_fourier, &
484  grid,leng1, table4, work, 0)
485 # endif
486 #else
487 # ifdef NAGFFT
488 ! NAG fft
489 ! will not allow float kind for nag
490  call error_handler ('fft_fourier_to_grid', &
491  'float kind not supported for nag fft')
492 
493  ! save input complex array in real format (herm.)
494  do k=1,lenc
495  data(:,k) = real(fourier(k,:))
496  enddo
497  do k=2,lenc-1
498  data(:,leng-k+2) = aimag(fourier(k,:))
499  enddo
500 
501  ! scale and transpose data
502  scale = sqrt(real(leng))
503  do j=1,num
504  grid(1:leng,j) = data(j,1:leng)*scale
505  enddo
506 # else
507 ! Temperton fft
508  do j=1,num
509  do k=1,lenc
510  data(2*k-1,j) = real(fourier(k,j))
511  data(2*k ,j) = aimag(fourier(k,j))
512  enddo
513  enddo
514  call fft991 (data,work,table99,ifax,1,leng2,leng,num,+1)
515  do j=1,num
516  grid(1:leng,j) = data(1:leng,j)
517  enddo
518 # endif
519 #endif
520 
521 !-----------------------------------------------------------------------
522 
523  end function fft_fourier_to_grid_float_2d
524 
525 
526 
527 !###############################################################################
528 !> \fn fft_grid_to_fourier_double_2d
529 !! \implements fft_grid_to_fourier
530 !!
531 !! \brief A function that returns multiple transforms in complex fourier space,
532 !! the first dimension must equal n/2+1 (where n is the size of each real
533 !! transform). The remaining dimensions must be the same size as the
534 !! input argument "grid".
535 !!
536 !! \code{.f90}
537 !! real (R8_KIND), intent(in), dimension(:,:) :: grid
538 !! complex(R8_KIND), dimension(lenc,size(grid,2)) :: fourier
539 !! \endcode
540 !!
541 !! \param [in] <grid> Multiple transforms in grid point space, the first
542 !! dimension must be n+1 (where n is the size of each real transform).
543 !!
544 !! \param [out] <fourier>
545 !!
546 !! \throw FATAL `"fft_init must be called"` \n
547 !! The initialization routine fft_init must be called before routines
548 !! fft_grid_to_fourier.
549 !!
550 !! \throw FATAL `"size of first dimension of input data is wrong"` \n
551 !! The real grid point field must have a first dimension equal to n+1
552 !! (where n is the size of each real transform). This message occurs
553 !! when using the SGI/Cray fft.
554 !!
555 !! \throw FATAL `"length of input data too small"` \n
556 !! The real grid point field must have a first dimension equal to n
557 !! (where n is the size of each real transform). This message occurs
558 !! when using the NAG or Temperton fft.
559 !!
560 !! \throw FATAL `"float kind not supported for nag fft"` \n
561 !! 32-bit real data is not supported when using the NAG fft. You
562 !! may try modifying this part of the code by uncommenting the
563 !! calls to the NAG library or less consider using the Temperton fft.
564 !!
565  function fft_grid_to_fourier_double_2d (grid) result (fourier)
567 !-----------------------------------------------------------------------
568 
569  real (R8_KIND), intent(in), dimension(:,:) :: grid
570  complex(R8_KIND), dimension(lenc,size(grid,2)) :: fourier
571 
572 #ifdef SGICRAY
573 # ifdef _CRAY
574 ! local storage for cray fft
575  real(R8_KIND), dimension((2*leng+4)*size(grid,2)) :: work
576 # else
577 ! local storage for sgi fft
578  real(R8_KIND), dimension(leng2) :: work
579 # endif
580 #else
581 # ifdef NAGFFT
582 ! local storage for nag fft
583  real(R8_KIND), dimension(size(grid,2),leng) :: data, work
584 # else
585 ! local storage for temperton fft
586  real, dimension(leng2,size(grid,2)) :: data
587  real, dimension(leng1,size(grid,2)) :: work
588 # endif
589 #endif
590 
591 #if defined(SGICRAY) || defined(NAGFFT)
592  real(R8_KIND) :: scale
593 #endif
594  integer :: j, k, num, len_grid
595 #ifdef NAGFFT
596  integer :: ifail
597 #endif
598 
599 !-----------------------------------------------------------------------
600 
601  if (.not.module_is_initialized) &
602  call error_handler ('fft_grid_to_fourier', &
603  'fft_init must be called.')
604 
605 !-----------------------------------------------------------------------
606 
607  len_grid = size(grid,1)
608 #ifdef SGICRAY
609  if (len_grid /= leng1) call error_handler ('fft_grid_to_fourier', &
610  'size of first dimension of input data is wrong')
611 #else
612  if (len_grid < leng) call error_handler ('fft_grid_to_fourier', &
613  'length of input data too small.')
614 #endif
615 !-----------------------------------------------------------------------
616 !----------------transform to fourier coefficients (+1)-----------------
617 
618  num = size(grid,2) ! number of transforms
619 #ifdef SGICRAY
620 ! Cray/SGI fft
621  scale = 1./float(leng)
622 # ifdef _CRAY
623  call scfftm (-1,leng,num,scale, grid,leng1, fourier,lenc, &
624  table8, work, 0)
625 # else
626  call dzfftm (-1,leng,num,scale, grid,leng1, fourier,lenc, &
627  table8, work, 0)
628 # endif
629 #else
630 # ifdef NAGFFT
631 ! NAG fft
632  do j=1,size(grid,2)
633  data(j,1:leng) = grid(1:leng,j)
634  enddo
635  call c06fpf ( num, leng, data, 's', table8, work, ifail )
636  scale = 1./sqrt(float(leng))
637  data = data * scale
638  fourier(1,:) = cmplx( data(:,1), 0. )
639  do k=2,lenc-1
640  fourier(k,:) = cmplx( data(:,k), data(:,leng-k+2) )
641  enddo
642  fourier(lenc,:) = cmplx( data(:,lenc), 0. )
643 # else
644 ! Temperton fft
645  do j=1,num
646  data(1:leng,j) = grid(1:leng,j)
647  enddo
648  call fft991 (data,work,table99,ifax,1,leng2,leng,num,-1)
649  do j=1,size(grid,2)
650  do k=1,lenc
651  fourier(k,j) = cmplx( data(2*k-1,j), data(2*k,j) )
652  enddo
653  enddo
654 # endif
655 #endif
656 !-----------------------------------------------------------------------
657 
658  end function fft_grid_to_fourier_double_2d
659 
660 
661 
662 !###############################################################################
663 !> \fn fft_fourier_to_grid_double_2d
664 !! \implements fft_fourier_to_grid
665 !!
666 !! \brief A function that returns multiple transforms in grid point space, the
667 !! first dimension must be n+1 (where n is the size of each real
668 !! transform). The remaining dimensions must be the same size as the
669 !! input argument "fourier".
670 !!
671 !! \code{.f90}
672 !! complex(R8_KIND), intent(in), dimension(:,:) :: fourier
673 !! real (R8_KIND), dimension(leng1,size(fourier,2)) :: grid
674 !!
675 !! \endcode
676 !!
677 !! \param [in] <fourier> Multiple transforms in complex fourier space, the
678 !! first dimension must equal n/2+1 (where n is the size of each
679 !! real transform).
680 !!
681 !! \param [out] <grid>
682 !!
683 !! \throw FATAL `"fft_init must be called"` \n
684 !! The initialization routine fft_init must be called before routines
685 !! fft_fourier_to_grid.
686 !!
687 !! \throw FATAL `"size of first dimension of input data is wrong"` \n
688 !! The complex Fourier field must have a first dimension equal to
689 !! n/2+1 (where n is the size of each real transform). This message
690 !! occurs when using the SGI/Cray fft.
691 !!
692 !! \throw FATAL `"length of input data too small"` \n
693 !! The complex Fourier field must have a first dimension greater
694 !! than or equal to n/2+1 (where n is the size of each real
695 !! transform). This message occurs when using the NAG or Temperton fft.
696 !!
697 !! \throw FATAL `"float kind not supported for nag fft"` \n
698 !! 32-bit real data is not supported when using the NAG fft. You may try
699 !! modifying this part of the code by uncommenting the calls to the NAG
700 !! library or less consider using the Temperton fft.
701 !!
702  function fft_fourier_to_grid_double_2d (fourier) result (grid)
704 !-----------------------------------------------------------------------
705 
706  complex(R8_KIND), intent(in), dimension(:,:) :: fourier
707  real (R8_KIND), dimension(leng1,size(fourier,2)) :: grid
708 
709 #ifdef SGICRAY
710 # ifdef _CRAY
711 ! local storage for cray fft
712  real(R8_KIND), dimension((2*leng+4)*size(fourier,2)) :: work
713 # else
714 ! local storage for sgi fft
715  real(R8_KIND), dimension(leng2) :: work
716 # endif
717 #else
718 # ifdef NAGFFT
719 ! local storage for nag fft
720  real(R8_KIND), dimension(size(fourier,2),leng) :: data, work
721 # else
722 ! local storage for temperton fft
723  real, dimension(leng2,size(fourier,2)) :: data
724  real, dimension(leng1,size(fourier,2)) :: work
725 # endif
726 #endif
727 
728 #if defined(SGICRAY) || defined(NAGFFT)
729  real(R8_KIND) :: scale
730 #endif
731  integer :: j, k, num, len_fourier
732 #ifdef NAGFFT
733  integer :: ifail
734 #endif
735 
736 !-----------------------------------------------------------------------
737 
738  if (.not.module_is_initialized) &
739  call error_handler ('fft_grid_to_fourier', &
740  'fft_init must be called.')
741 
742 !-----------------------------------------------------------------------
743 
744  len_fourier = size(fourier,1)
745  num = size(fourier,2) ! number of transforms
746 
747 #ifdef SGICRAY
748  if (len_fourier /= lenc) call error_handler ('fft_fourier_to_grid', &
749  'size of first dimension of input data is wrong')
750 #else
751  if (len_fourier < lenc) call error_handler ('fft_fourier_to_grid', &
752  'length of input data too small.')
753 #endif
754 !-----------------------------------------------------------------------
755 !----------------inverse transform to real space (-1)-------------------
756 
757 #ifdef SGICRAY
758 ! Cray/SGI fft
759  scale = 1.0
760 # ifdef _CRAY
761  call csfftm (+1,leng,num,scale, fourier,len_fourier, &
762  grid,leng1, table8, work, 0)
763 # else
764  call zdfftm (+1,leng,num,scale, fourier,len_fourier, &
765  grid,leng1, table8, work, 0)
766 # endif
767 #else
768 # ifdef NAGFFT
769 ! NAG fft
770 
771  ! save input complex array in real format (herm.)
772  do k=1,lenc
773  data(:,k) = real(fourier(k,:))
774  enddo
775  do k=2,lenc-1
776  data(:,leng-k+2) = aimag(fourier(k,:))
777  enddo
778 
779  call c06gqf ( num, leng, data, ifail )
780  call c06fqf ( num, leng, data, 's', table8, work, ifail )
781 
782  ! scale and transpose data
783  scale = sqrt(real(leng))
784  do j=1,num
785  grid(1:leng,j) = data(j,1:leng)*scale
786  enddo
787 # else
788 ! Temperton fft
789  do j=1,num
790  do k=1,lenc
791  data(2*k-1,j) = real(fourier(k,j))
792  data(2*k ,j) = aimag(fourier(k,j))
793  enddo
794  enddo
795  call fft991 (data,work,table99,ifax,1,leng2,leng,num,+1)
796  do j=1,num
797  grid(1:leng,j) = data(1:leng,j)
798  enddo
799 # endif
800 #endif
801 
802 !-----------------------------------------------------------------------
803 
804  end function fft_fourier_to_grid_double_2d
805 
806 !###############################################################################
807 ! interface overloads
808 !###############################################################################
809 
810 
811 
812 !###############################################################################
813 !> \fn fft_grid_to_fourier_float_3d
814 !! \implements fft_grid_to_fourier
815 !!
816 !! \code{.f90}
817 !! real (R4_KIND), intent(in), dimension(:,:,:) :: grid
818 !! complex(R4_KIND), dimension(lenc,size(grid,2),size(grid,3)) :: fourier
819 !! \endcode
820 !!
821 !! \param [in] <grid>
822 !!
823 !! \param [out] <fourier>
824 !!
825  function fft_grid_to_fourier_float_3d (grid) result (fourier)
827 !-----------------------------------------------------------------------
828  real (R4_KIND), intent(in), dimension(:,:,:) :: grid
829  complex(R4_KIND), dimension(lenc,size(grid,2),size(grid,3)) :: fourier
830  integer :: n
831 !-----------------------------------------------------------------------
832 
833  do n = 1, size(grid,3)
834  fourier(:,:,n) = fft_grid_to_fourier_float_2d(grid(:,:,n))
835  enddo
836 
837 !-----------------------------------------------------------------------
838 
839  end function fft_grid_to_fourier_float_3d
840 
841 
842 
843 !###############################################################################
844 !> \fn fft_fourier_to_grid_float_3d
845 !! \implements fft_fourier_to_grid
846 !!
847 !! \code{.f90}
848 !! complex(R4_KIND), intent(in), dimension(:,:,:) :: fourier
849 !! real (R4_KIND), dimension(leng1,size(fourier,2),size(fourier,3)) :: grid
850 !! \endcode
851 !!
852 !! \param [in] <fourier>
853 !!
854 !! \param [out] <grid>
855 !!
856  function fft_fourier_to_grid_float_3d (fourier) result (grid)
858 !-----------------------------------------------------------------------
859  complex(R4_KIND), intent(in), dimension(:,:,:) :: fourier
860  real (R4_KIND), dimension(leng1,size(fourier,2),size(fourier,3)) :: grid
861  integer :: n
862 !-----------------------------------------------------------------------
863 
864  do n = 1, size(fourier,3)
865  grid(:,:,n) = fft_fourier_to_grid_float_2d(fourier(:,:,n))
866  enddo
867 
868 !-----------------------------------------------------------------------
869 
870  end function fft_fourier_to_grid_float_3d
871 
872 
873 
874 !###############################################################################
875 !> \fn fft_grid_to_fourier_double_3d
876 !! \implements fft_grid_to_fourier
877 !!
878 !! \code{.f90}
879 !! real (R8_KIND), intent(in), dimension(:,:,:) :: grid
880 !! complex(R8_KIND), dimension(lenc,size(grid,2),size(grid,3)) :: fourier
881 !! \endcode
882 !!
883 !! \param [in] <grid>
884 !!
885 !! \param [out] <fourier>
886 !!
887  function fft_grid_to_fourier_double_3d (grid) result (fourier)
889 !-----------------------------------------------------------------------
890  real (R8_KIND), intent(in), dimension(:,:,:) :: grid
891  complex(R8_KIND), dimension(lenc,size(grid,2),size(grid,3)) :: fourier
892  integer :: n
893 !-----------------------------------------------------------------------
894 
895  do n = 1, size(grid,3)
896  fourier(:,:,n) = fft_grid_to_fourier_double_2d(grid(:,:,n))
897  enddo
898 
899 !-----------------------------------------------------------------------
900 
901  end function fft_grid_to_fourier_double_3d
902 
903 
904 
905 !#######################################################################
906 !> \fn fft_fourier_to_grid_double_3d
907 !! \implements fft_fourier_to_grid
908 !!
909 !! \param [in] <fourier>
910 !!
911 !! \param [out] <grid>
912 !!
913 !! \code{.f90}
914 !! complex(R8_KIND), intent(in), dimension(:,:,:) :: fourier
915 !!
916 !! real(R8_KIND), dimension(leng1,size(fourier,2),size(fourier,3)) :: grid
917 !! \endcode
918 !!
919 !! \note The original documentation for this method contradicts the actual code.
920 !! The type of `fourier` is `complex(R8_KIND)` while the type of `grid` is
921 !! `real(R8_KIND)`. Here is the original documentation:
922 !! \code
923 !! <FUNCTION NAME="fft_fourier_to_grid_double_3d" INTERFACE="fft_fourier_to_grid">
924 !! <IN NAME="fourier" TYPE="real(R8_KIND)" DIM="(:,:,:)"></IN>
925 !! <OUT NAME="grid" TYPE="complex(R8_KIND)" DIM="(leng1,size(fourier,2),size(fourier,3))"> </OUT>
926 !! </FUNCTION>
927 !! \endcode
928 !!
929  function fft_fourier_to_grid_double_3d (fourier) result (grid)
931 !-----------------------------------------------------------------------
932  complex(R8_KIND), intent(in), dimension(:,:,:) :: fourier
933  real (R8_KIND), dimension(leng1,size(fourier,2),size(fourier,3)) :: grid
934  integer :: n
935 !-----------------------------------------------------------------------
936 
937  do n = 1, size(fourier,3)
938  grid(:,:,n) = fft_fourier_to_grid_double_2d(fourier(:,:,n))
939  enddo
940 
941 !-----------------------------------------------------------------------
942 
943  end function fft_fourier_to_grid_double_3d
944 
945 
946 
947 !###############################################################################
948 !> \fn fft_init
949 !! \ingroup fft_mod
950 !!
951 !! \brief This routine must be called to initialize the size of a
952 !! single transform and setup trigonometric constants.
953 !!
954 !! This routine must be called once to initialize the size of a
955 !! single transform. To change the size of the transform the
956 !! routine fft_exit must be called before re-initialing with fft_init.
957 !!
958 !!
959 !! \code{.f90}
960 !! integer, intent(in) :: n
961 !! \endcode
962 !!
963 !! \param [in] <n> The number of real values in a single sequence of data.
964 !! The resulting transformed data will have n/2+1 pairs of
965 !! complex values.
966 !!
967 !! \throw FATAL `"attempted to reinitialize fft"` \n
968 !! You must call fft_exit before calling fft_init for a second time.
969 !!
970 !! <TEMPLATE>
971 !! call fft_init ( n )
972 !! </TEMPLATE>
973  subroutine fft_init (n)
975 !-----------------------------------------------------------------------
976  integer, intent(in) :: n
977 
978 #ifdef SGICRAY
979  real (R4_KIND) :: dummy4(1)
980  complex(R4_KIND) :: cdummy4(1)
981  real (R8_KIND) :: dummy8(1)
982  complex(R8_KIND) :: cdummy8(1)
983  integer :: isys(0:1)
984 #else
985 # ifdef NAGFFT
986  real(R8_KIND) :: data8(n), work8(n)
987  real(R4_KIND) :: data4(n), work4(n)
988  integer :: ifail4, ifail8
989 # endif
990 #endif
991 !-----------------------------------------------------------------------
992 ! --- fourier transform initialization ----
993 
994  if (module_is_initialized) &
995  call error_handler ('fft_init', 'attempted to reinitialize fft')
996 
997 ! write file version to log file
998  if (do_log) then
999  call write_version_number("FFT_MOD", version)
1000  do_log = .false.
1001  endif
1002 
1003 ! variables that save length of transform
1004  leng = n; leng1 = n+1; leng2 = n+2; lenc = n/2+1
1005 
1006 #ifdef SGICRAY
1007 # ifdef _CRAY
1008 ! initialization for cray
1009 ! float kind may not apply for cray
1010  allocate (table4(100+2*leng), table8(100+2*leng)) ! size may be too large?
1011  call scfftm (0,leng,1,0.0, dummy4, 1, cdummy4, 1, table4, dummy4, 0)
1012  call scfftm (0,leng,1,0.0, dummy8, 1, cdummy8, 1, table8, dummy8, 0)
1013 # else
1014 ! initialization for sgi
1015  allocate (table4(leng+256), table8(leng+256))
1016  isys(0) = 1
1017  call scfftm (0,leng,1,0.0, dummy4, 1, cdummy4, 1, table4, dummy8, isys)
1018  call dzfftm (0,leng,1,0.0, dummy8, 1, cdummy8, 1, table8, dummy8, isys)
1019 # endif
1020 #else
1021 # ifdef NAGFFT
1022 ! initialization for nag fft
1023  ifail8 = 0
1024  allocate (table8(100+2*leng)) ! size may be too large?
1025  call c06fpf ( 1, leng, data8, 'i', table8, work8, ifail8 )
1026 
1027 ! will not allow float kind for nag
1028  ifail4 = 0
1029 !!!!! allocate (table4(100+2*leng))
1030 !!!!! call c06fpe ( 1, leng, data4, 'i', table4, work4, ifail4 )
1031 
1032  if (ifail4 /= 0 .or. ifail8 /= 0) then
1033  call error_handler ('fft_init', 'nag fft initialization error')
1034  endif
1035 # else
1036 ! initialization for Temperton fft
1037  allocate (table99(3*leng/2+1))
1038  allocate (ifax(10))
1039  call set99 ( table99, ifax, leng )
1040 # endif
1041 #endif
1042 
1043  module_is_initialized = .true.
1044 
1045 !-----------------------------------------------------------------------
1046 
1047  end subroutine fft_init
1048 
1049 
1050 
1051 !###############################################################################
1052 !> \fn fft_end
1053 !! \ingroup fft_mod
1054 !!
1055 !! \brief This routine is called to unset the transform size and deallocate
1056 !! memory.
1057 !!
1058 !! This routine is called to unset the transform size and
1059 !! deallocate memory. It can not be called unless fft_init
1060 !! has already been called. There are no arguments.
1061 !!
1062 !! \throw FATAL `"attempt to un-initialize fft that has not been initialized"` \n
1063 !! You can not call fft_end unless fft_init has been called.
1064 !!
1065 !! <TEMPLATE>
1066 !! call fft_end
1067 !! </TEMPLATE>
1068  subroutine fft_end
1070 ! --- fourier transform un-initialization ----
1071 
1072  if (.not.module_is_initialized) &
1073  call error_handler ('fft_end', &
1074  'attempt to un-initialize fft that has not been initialized')
1075 
1076  leng = 0; leng1 = 0; leng2 = 0; lenc = 0
1077 
1078  if (allocated(table4)) deallocate (table4)
1079  if (allocated(table8)) deallocate (table8)
1080  if (allocated(table99)) deallocate (table99)
1081 
1082  module_is_initialized = .false.
1083 
1084 !-----------------------------------------------------------------------
1085 
1086  end subroutine fft_end
1087 
1088 
1089 
1090 !###############################################################################
1091 !> \fn error_handler
1092 !! \ingroup fft_mod
1093 !!
1094 !! \brief Wrapper for handling errors
1095 !!
1096  subroutine error_handler ( routine, message )
1097  character(len=*), intent(in) :: routine, message
1098 
1099  call error_mesg ( routine, message, fatal )
1100 
1101 ! print *, 'ERROR: ',trim(routine)
1102 ! print *, 'ERROR: ',trim(message)
1103 ! stop 111
1104 
1105  end subroutine error_handler
1106 
1107 
1108 
1109 end module fft_mod
1110 
1111 #ifdef test_fft
1112 program test
1113 use fft_mod
1114 integer, parameter :: lot = 2
1115 real , allocatable :: ain(:,:), aout(:,:)
1116 complex, allocatable :: four(:,:)
1117 integer :: i, j, m, n
1118 integer :: ntrans(2) = (/ 60, 90 /)
1119 
1120 ! test multiple transform lengths
1121  do m = 1,2
1122 
1123  ! set up input data
1124  n = ntrans(m)
1125  allocate (ain(n+1,lot),aout(n+1,lot),four(n/2+1,lot))
1126  call random_number (ain(1:n,:))
1127  aout(1:n,:) = ain(1:n,:)
1128 
1129  call fft_init (n)
1130  ! transform grid to fourier and back
1131  four = fft_grid_to_fourier(aout)
1132  aout = fft_fourier_to_grid(four)
1133 
1134  ! print original and transformed
1135  do j=1,lot
1136  do i=1,n
1137  write (*,'(2i4,3(2x,f15.9))') j, i, ain(i,j), aout(i,j), aout(i,j)-ain(i,j)
1138  enddo
1139  enddo
1140 
1141  call fft_end
1142  deallocate (ain,aout,four)
1143  enddo
1144 
1145 end program test
1146 #endif
1147 
1148 !> @} end of the fft_mod module
1149 
1150 ! <INFO>
1151 ! <REFERENCE>
1152 ! For the SGI/Cray version refer to the manual pages for
1153 ! DZFFTM, ZDFFTM, SCFFTM, and CSFFTM.
1154 ! </REFERENCE>
1155 ! <REFERENCE>
1156 ! For the NAG version refer to the NAG documentation for
1157 ! routines C06FPF, C06FQF, and C06GQF.
1158 ! </REFERENCE>
1159 ! <PRECOMP FLAG="-D NAGFFT">
1160 ! -D NAGFFT
1161 ! On non-Cray/SGI machines, set to use the NAG library FFT routines.
1162 ! Otherwise the Temperton FFT is used by default.
1163 ! </PRECOMP>
1164 ! <PRECOMP FLAG="-D test_fft">
1165 ! Provides source code for a simple test program.
1166 ! The program generates several sequences of real data.
1167 ! This data is transformed to Fourier space and back to real data,
1168 ! then compared to the original real data.
1169 ! </PRECOMP>
1170 ! <LOADER FLAG="-lscs">
1171 ! On SGI machines the scientific library needs to be loaded by
1172 ! linking with:
1173 ! </LOADER>
1174 ! <LOADER FLAG="-L/usr/local/lib -lnag">
1175 ! If using the NAG library, the following loader options (or
1176 ! something similar) may be necessary:
1177 ! </LOADER>
1178 ! <NOTE>
1179 ! The routines are overloaded for 2d and 3d versions.
1180 ! The 2d versions copy data into 3d arrays then calls the 3d interface.
1181 !
1182 ! On SGI/Cray machines:
1183 !
1184 ! There are single (32-bit) and full (64-bit) versions.
1185 ! For Cray machines the single precision version does not apply.
1186 !
1187 ! On non-SGI/CRAY machines:
1188 !
1189 ! The NAG library option uses the "full" precision NAG
1190 ! routines (C06FPF,C06FQF,C06GQF). Users may have to specify
1191 ! a 64-bit real compiler option (e.g., -r8).
1192 !
1193 ! The stand-alone Temperton FFT option works for the
1194 ! real precision specified at compile time.
1195 ! If you compiled with single (32-bit) real precision
1196 ! then FFT's cannot be computed at full (64-bit) precision.
1197 ! </NOTE>
1198 ! </INFO>
Definition: fms.F90:20
real(r4_kind) function, dimension(leng1, size(fourier, 2), size(fourier, 3)) fft_fourier_to_grid_float_3d(fourier)
Definition: fft.F90:857
subroutine, public fft991(a, work, trigs, ifax, inc, jump, n, lot, isign)
Definition: fft99.F90:578
subroutine error_handler(routine, message)
Definition: fft.F90:1097
real(r8_kind) function, dimension(leng1, size(fourier, 2), size(fourier, 3)) fft_fourier_to_grid_double_3d(fourier)
Definition: fft.F90:930
integer, parameter r4_kind
Definition: platform.F90:24
logical module_is_initialized
Definition: fft.F90:240
real(r4_kind) function, dimension(leng1, size(fourier, 2)) fft_fourier_to_grid_float_2d(fourier)
Definition: fft.F90:425
integer, dimension(:), allocatable ifax
Definition: fft.F90:235
complex(r4_kind) function, dimension(lenc, size(grid, 2)) fft_grid_to_fourier_float_2d(grid)
Definition: fft.F90:289
subroutine, public fft_init(n)
Definition: fft.F90:974
integer lenc
Definition: fft.F90:238
subroutine, public set99(trigs, ifax, n)
Definition: fft99.F90:731
Fortran module for Eigen FFTs.
Definition: fft.F90:24
integer leng2
Definition: fft.F90:238
complex(r8_kind) function, dimension(lenc, size(grid, 2)) fft_grid_to_fourier_double_2d(grid)
Definition: fft.F90:566
real(r8_kind), dimension(:), allocatable table8
Definition: fft.F90:232
integer, parameter r8_kind
Definition: platform.F90:24
integer leng1
Definition: fft.F90:238
real(r4_kind), dimension(:), allocatable table4
Definition: fft.F90:233
complex(r8_kind) function, dimension(lenc, size(grid, 2), size(grid, 3)) fft_grid_to_fourier_double_3d(grid)
Definition: fft.F90:888
real, dimension(:), allocatable table99
Definition: fft.F90:234
complex(r4_kind) function, dimension(lenc, size(grid, 2), size(grid, 3)) fft_grid_to_fourier_float_3d(grid)
Definition: fft.F90:826
logical do_log
Definition: fft.F90:237
subroutine, public fft_end
This routine is called to unset the transform size and deallocate memory.
Definition: fft.F90:1069
subroutine, public error_mesg(routine, message, level)
Definition: fms.F90:529
real(r8_kind) function, dimension(leng1, size(fourier, 2)) fft_fourier_to_grid_double_2d(fourier)
Definition: fft.F90:703
integer leng
Definition: fft.F90:238