tab_calc_fftw3_real_bwd Subroutine

public subroutine tab_calc_fftw3_real_bwd(tab_in, tab_ou, long, larg)

Note

Subroutine that transforms backward a real 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
complex(kind=R8), intent(in), dimension(1:long, 1:larg) :: tab_in

array to transform

real(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_real_bwd~~CallsGraph proc~tab_calc_fftw3_real_bwd tab_calc_fftw3_real_bwd interface~fftw_execute_dft_c2r fftw_execute_dft_c2r proc~tab_calc_fftw3_real_bwd->interface~fftw_execute_dft_c2r omp_get_thread_num omp_get_thread_num proc~tab_calc_fftw3_real_bwd->omp_get_thread_num

Called by

proc~~tab_calc_fftw3_real_bwd~~CalledByGraph proc~tab_calc_fftw3_real_bwd tab_calc_fftw3_real_bwd program~test_fftw3 test_fftw3 program~test_fftw3->proc~tab_calc_fftw3_real_bwd

Source Code

   subroutine tab_calc_fftw3_real_bwd(tab_in, tab_ou, long, larg)
   implicit none
   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*
   real   (kind=R8), dimension(1:long, 1:larg), intent(out) :: tab_ou   !! *transformed array*

      integer(kind=I4) :: ithread

      ithread = omp_get_thread_num()

      tab_cmp_b_i(ithread)%tab(1:long, 1:larg) = tab_in(1:long, 1:larg)

      call fftw_execute_dft_c2r(tab_plan_b(ithread), tab_cmp_b_i(ithread)%tab(1:long, 1:larg),   &  !
                                                     tab_rea_b_o(ithread)%tab(1:long, 1:larg))      !

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

   return
   endsubroutine tab_calc_fftw3_real_bwd