calc_fftw3_real_bwd Subroutine

public subroutine calc_fftw3_real_bwd(tab_in, tab_ou, long, larg, planner_flag)

Note

Subroutine that transforms backward a double real array. For speed reasons FFTW will always work on the same memory area, until the plans are destroyed of course. 1 FFT distributed on several threads

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

integer(kind=I4), intent(in), optional :: planner_flag

planning option, FFTW_ESTIMATE for example


Calls

proc~~calc_fftw3_real_bwd~~CallsGraph proc~calc_fftw3_real_bwd calc_fftw3_real_bwd interface~fftw_execute_dft_c2r fftw_execute_dft_c2r proc~calc_fftw3_real_bwd->interface~fftw_execute_dft_c2r interface~fftw_plan_with_nthreads fftw_plan_with_nthreads proc~calc_fftw3_real_bwd->interface~fftw_plan_with_nthreads omp_get_num_procs omp_get_num_procs proc~calc_fftw3_real_bwd->omp_get_num_procs proc~end_fftw3 end_fftw3 proc~calc_fftw3_real_bwd->proc~end_fftw3 proc~init_fftw3_real init_fftw3_real proc~calc_fftw3_real_bwd->proc~init_fftw3_real proc~desalloc_fftw3 desalloc_fftw3 proc~end_fftw3->proc~desalloc_fftw3 proc~destroy_plan_fftw3 destroy_plan_fftw3 proc~end_fftw3->proc~destroy_plan_fftw3 proc~alloc_fftw3_real alloc_fftw3_real proc~init_fftw3_real->proc~alloc_fftw3_real proc~make_plan_fftw3_real make_plan_fftw3_real proc~init_fftw3_real->proc~make_plan_fftw3_real interface~fftw_alloc_complex fftw_alloc_complex proc~alloc_fftw3_real->interface~fftw_alloc_complex interface~fftw_alloc_real fftw_alloc_real proc~alloc_fftw3_real->interface~fftw_alloc_real interface~fftw_free fftw_free proc~desalloc_fftw3->interface~fftw_free interface~fftw_destroy_plan fftw_destroy_plan proc~destroy_plan_fftw3->interface~fftw_destroy_plan interface~fftw_plan_dft_c2r_2d fftw_plan_dft_c2r_2d proc~make_plan_fftw3_real->interface~fftw_plan_dft_c2r_2d interface~fftw_plan_dft_r2c_2d fftw_plan_dft_r2c_2d proc~make_plan_fftw3_real->interface~fftw_plan_dft_r2c_2d

Called by

proc~~calc_fftw3_real_bwd~~CalledByGraph proc~calc_fftw3_real_bwd calc_fftw3_real_bwd program~test_fftw3 test_fftw3 program~test_fftw3->proc~calc_fftw3_real_bwd

Source Code

   subroutine calc_fftw3_real_bwd(tab_in, tab_ou, long, larg, planner_flag)
   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), intent(in),                optional    :: planner_flag   !! *planning option, [[fftw3(module):FFTW_ESTIMATE]] for example*

      integer(kind=I4) :: plan_flag

      if ( .not. present(planner_flag) ) then

         plan_flag = FFTW_ESTIMATE

      else

         plan_flag = planner_flag

      endif

      if ( any( FFT_DIM(1:2) /= [long, larg] ) ) then

         if ( sum(FFT_DIM(1:2)) /= 0 ) call end_fftw3()

         call fftw_plan_with_nthreads(nthreads = omp_get_num_procs())
         call init_fftw3_real(long = long, larg = larg, plan_flag = plan_flag)

      endif

      cmp_b_i(1:long, 1:larg) = tab_in(1:long, 1:larg)

      call fftw_execute_dft_c2r(plan_b, cmp_b_i(1:long, 1:larg), rea_b_o(1:long, 1:larg))

      tab_ou(1:long, 1:larg) = rea_b_o(1:long, 1:larg) / sqrt(real(long*larg, kind=r8))

   return
   endsubroutine calc_fftw3_real_bwd