tab_calc_fftw3 Subroutine

public subroutine tab_calc_fftw3(sens, tab_in, tab_ou, long, larg)

Note

Subroutine that transforms forward or bacward a double complex array. For speed reasons FFTW will always work on the same memory area, until the plans are destroyed of course. several FFT on single thread each

Arguments

Type IntentOptional Attributes Name
integer(kind=I4), intent(in) :: sens

=FORWARD or =BACKWARD

complex(kind=R8), intent(in), dimension(1:long, 1:larg) :: tab_in

array to transform

complex(kind=R8), intent(out), dimension(1:long, 1:larg) :: tab_ou

transformed array

integer(kind=I4), intent(in) :: long

first 2D array dimension

integer(kind=I4), intent(in) :: larg

second 2D array dimension


Calls

proc~~tab_calc_fftw3~~CallsGraph proc~tab_calc_fftw3 tab_calc_fftw3 interface~fftw_execute_dft fftw_execute_dft proc~tab_calc_fftw3->interface~fftw_execute_dft omp_get_thread_num omp_get_thread_num proc~tab_calc_fftw3->omp_get_thread_num

Called by

proc~~tab_calc_fftw3~~CalledByGraph proc~tab_calc_fftw3 tab_calc_fftw3 program~test_fftw3 test_fftw3 program~test_fftw3->proc~tab_calc_fftw3

Source Code

   subroutine tab_calc_fftw3(sens, tab_in, tab_ou, long, larg)
   implicit none
   integer(kind=I4), intent(in ) :: sens                                !! *```=FORWARD``` or ```=BACKWARD```*
   integer(kind=I4), intent(in ) :: long                                !! *first  2D array dimension*
   integer(kind=I4), intent(in ) :: larg                                !! *second 2D array dimension*
   complex(kind=R8), dimension(1:long, 1:larg), intent(in ) :: tab_in   !! *array to transform*
   complex(kind=R8), dimension(1:long, 1:larg), intent(out) :: tab_ou   !! *transformed array*

      integer(kind=I4) :: ithread

      ithread = omp_get_thread_num()
      select case(sens)

         case(FORWARD)
            tab_cmp_f_i(ithread)%tab(1:long, 1:larg) = tab_in(1:long, 1:larg)
            call fftw_execute_dft(tab_plan_f(ithread), tab_cmp_f_i(ithread)%tab(1:long, 1:larg),   &  !
                                                       tab_cmp_f_o(ithread)%tab(1:long, 1:larg))      !

            tab_ou(1:long, 1:larg) = tab_cmp_f_o(ithread)%tab(1:long, 1:larg)/sqrt(real(long*larg, kind=r8))

         case(BACKWARD)
            tab_cmp_b_i(ithread)%tab(1:long, 1:larg) = tab_in(1:long, 1:larg)
            call fftw_execute_dft(tab_plan_b(ithread), tab_cmp_b_i(ithread)%tab(1:long, 1:larg),   &  !
                                                       tab_cmp_b_o(ithread)%tab(1:long, 1:larg))      !

            tab_ou(1:long, 1:larg) = tab_cmp_b_o(ithread)%tab(1:long, 1:larg)/sqrt(real(long*larg, kind=r8))

      endselect

   return
   endsubroutine tab_calc_fftw3