Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Derived types and interfaces | Functions/Subroutines | Variables
sll_m_fft Module Reference

Description

FFT interface for FFTW.

Author
Katharina Kormann (IPP)
Yaman Güçlü (IPP)

Complex-to-Complex DFT (c2c)

The complex coefficients returned by the 'c2c' discrete Fourier transform correspond to the interval \([0,2\pi)\), which is not symmetric with respect to 0; nevertheless, because of the \(2\pi\)-periodicity of the transform, the coefficients in \((\pi,2\pi)\) are identical to those in \((-\pi,0)\). Therefore, if we think in terms of positive and negative frequencies, the \(N\) complex coefficients returned by the 'c2c' transform are ordered as

\[ \left[ c_0, c_1, c_2,..., c_{N/2}, c_{(-N+N/2+1)}, ..., c_{-2}, c_{-1} \right], \]

where the integer division is rounded down toward zero, so that

\[ -N+N/2+1 = \begin{cases} -N/2+1 & \text{if $N$ is even}, \\ -N/2 & \text{if $N$ is odd}. \end{cases} \]

Real-to-Complex DFT (r2c)

If the input data to be transformed is real, then we know that the Fourier transform has real and imaginary parts that are respectively symmetric and anti-symmetric around 0; therefore all the necessary information is contained in the interval \([0,\pi]\). The 'r2c' discrete Fourier transform takes advantage of this by computing only the \(N/2+1\) complex Fourier coefficients

\[ \left[ c_0, c_1, c_2, ..., c_{N/2} \right], \]

where \(c_0\) has zero imaginary part, as well as \(c_{N/2}\) if N is even.

Real-to-Real DFT (r2r)

When the input data array is real it is sometimes more convenient to store the complex Fourier coefficients in a real output array of identical size. The 'r2r' discrete Fourier transform provides this ability by using the so-called 'halcomplex' format, where the \(N/2+1\) complex Fourier coefficients are given in the form of \(N\) real values

\[ \left[ r_0, r_1, r_2, ..., r_{N/2}, i_{(N+1)/2-1}, ..., i_{-2}, i_{-1} \right]. \]

Here \(r_k\) and \(i_k\) are respectively the real and imaginary parts of the complex coefficient \(c_k\).

Derived types and interfaces

type  sll_t_fft
 Type for FFT plan in SeLaLib. More...
 

Functions/Subroutines

subroutine, public sll_s_fft_print_defaultlib ()
 Function to print the FFT library behind the interface. More...
 
complex(c_double_complex) function, dimension(:), pointer, public sll_f_fft_allocate_aligned_complex (n)
 Function to allocate an aligned complex array. More...
 
real(c_double) function, dimension(:), pointer, public sll_f_fft_allocate_aligned_real (n)
 Function to allocate an aligned real array. More...
 
subroutine, public sll_s_fft_deallocate_aligned_complex (data)
 Function to deallocate an aligned complex array. More...
 
subroutine, public sll_s_fft_deallocate_aligned_real (data)
 Function to deallocate an aligned real array. More...
 
subroutine, public sll_s_fft_get_k_list_c2c_1d (plan, k_list)
 Get a list of the k modes in a c2c FFT k_list = [0, 1, 2, ..., n/2, (n/2+1)-n, ..., -2, -1]. More...
 
subroutine, public sll_s_fft_get_k_list_r2r_1d (plan, k_list)
 Get a list of the k modes in an r2r FFT k_list = [0, 1, 2, ..., n/2, (n+1)/2-1, ..., 2, 1]. More...
 
subroutine, public sll_s_fft_get_k_list_r2c_1d (plan, k_list)
 Get a list of the k modes in an r2c FFT k_list = [0, 1, 2, ..., n/2]. More...
 
complex(kind=f64) function, public sll_f_fft_get_mode_r2c_1d (plan, data, k)
 Function to reconstruct a complex FFT mode from the data of a r2r transform. More...
 
subroutine, public sll_s_fft_set_mode_c2r_1d (plan, data, ck, k)
 Subroutine to set a complex mode to the real representation of r2r. More...
 
subroutine, public sll_s_fft_init_c2c_1d (plan, nx, array_in, array_out, direction, normalized, aligned, optimization)
 Create new 1d complex to complex plan. More...
 
subroutine, public sll_s_fft_exec_c2c_1d (plan, array_in, array_out)
 Compute fast Fourier transform in complex to complex mode. More...
 
subroutine, public sll_s_fft_init_c2c_2d (plan, nx, ny, array_in, array_out, direction, normalized, aligned, optimization)
 Create new 2d complex to complex plan. More...
 
subroutine, public sll_s_fft_exec_c2c_2d (plan, array_in, array_out)
 Compute fast Fourier transform in complex to complex mode. More...
 
subroutine, public sll_s_fft_init_c2c_3d (plan, nx, ny, nz, array_in, array_out, direction, normalized, aligned, optimization)
 Create new 3d complex to complex plan. More...
 
subroutine, public sll_s_fft_exec_c2c_3d (plan, array_in, array_out)
 Compute fast Fourier transform in complex to complex mode. More...
 
subroutine, public sll_s_fft_init_r2r_1d (plan, nx, array_in, array_out, direction, normalized, aligned, optimization)
 Create new 1d real to real plan. More...
 
subroutine, public sll_s_fft_exec_r2r_1d (plan, array_in, array_out)
 Compute fast Fourier transform in real to real mode. More...
 
subroutine, public sll_s_fft_init_r2c_1d (plan, nx, array_in, array_out, normalized, aligned, optimization)
 Create new 1d real to complex plan for forward FFT. More...
 
subroutine, public sll_s_fft_exec_r2c_1d (plan, array_in, array_out)
 Compute fast Fourier transform in real to complex mode. More...
 
subroutine, public sll_s_fft_init_r2c_2d (plan, nx, ny, array_in, array_out, normalized, aligned, optimization)
 Create new 2d complex to real plan for forward FFT. More...
 
subroutine, public sll_s_fft_exec_r2c_2d (plan, array_in, array_out)
 Compute fast Fourier transform in real to complex mode. More...
 
subroutine, public sll_s_fft_init_c2r_1d (plan, nx, array_in, array_out, normalized, aligned, optimization)
 Create new 1d complex to real plan for backward FFT. More...
 
subroutine, public sll_s_fft_exec_c2r_1d (plan, array_in, array_out)
 Compute fast Fourier transform in complex to real mode. More...
 
subroutine, public sll_s_fft_init_c2r_2d (plan, nx, ny, array_in, array_out, normalized, aligned, optimization)
 Create new 2d real to complex plan for backward FFT. More...
 
subroutine, public sll_s_fft_exec_c2r_2d (plan, array_in, array_out)
 Compute fast Fourier transform in complex to real mode. More...
 
subroutine, public sll_s_fft_free (plan)
 Delete a plan. More...
 

Variables

integer, parameter, public sll_p_fft_forward = FFTW_FORWARD
 Flag to specify forward FFT (negative sign) [value -1]. More...
 
integer, parameter, public sll_p_fft_backward = FFTW_BACKWARD
 Flag to specify backward FFT (positive sign) [value 1]. More...
 
integer(c_int), parameter, public sll_p_fft_measure = FFTW_MEASURE
 FFTW planning-rigor flag FFTW_MEASURE (optimized plan) NOTE: planner overwrites the input array during planning [value 0]. More...
 
integer(c_int), parameter, public sll_p_fft_patient = FFTW_PATIENT
 FFTW planning-rigor flag FFTW_PATIENT (more optimizaton than MEASURE) NOTE: planner overwrites the input array during planning [value 32]. More...
 
integer(c_int), parameter, public sll_p_fft_estimate = FFTW_ESTIMATE
 FFTW planning-rigor flag FFTW_ESTIMATE (simple heuristic for planer) [value 64]. This is our default value. More...
 
integer(c_int), parameter, public sll_p_fft_exhaustive = FFTW_EXHAUSTIVE
 FFTW planning-rigor flag FFTW_EXHAUSTIVE (more optimization than PATIENT) NOTE: planner overwrites the input array during planning [value 8]. More...
 
integer(c_int), parameter, public sll_p_fft_wisdom_only = FFTW_WISDOM_ONLY
 
integer, parameter p_fftw_c2c_1d = 0
 
integer, parameter p_fftw_r2r_1d = 1
 
integer, parameter p_fftw_r2c_1d = 2
 
integer, parameter p_fftw_c2r_1d = 3
 
integer, parameter p_fftw_c2c_2d = 4
 
integer, parameter p_fftw_r2r_2d = 5
 
integer, parameter p_fftw_r2c_2d = 6
 
integer, parameter p_fftw_c2r_2d = 7
 
integer, parameter p_fftw_c2c_3d = 8
 

Function/Subroutine Documentation

◆ sll_f_fft_allocate_aligned_complex()

complex(c_double_complex) function, dimension(:), pointer, public sll_m_fft::sll_f_fft_allocate_aligned_complex ( integer(kind=i32), intent(in)  n)

Function to allocate an aligned complex array.

Parameters
[in]nSize of the pointer
Returns
Array to be allocated

Definition at line 177 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_f_fft_allocate_aligned_real()

real(c_double) function, dimension(:), pointer, public sll_m_fft::sll_f_fft_allocate_aligned_real ( integer(kind=i32), intent(in)  n)

Function to allocate an aligned real array.

Parameters
[in]nSize of the pointer
Returns
Array to be allocated

Definition at line 196 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_f_fft_get_mode_r2c_1d()

complex(kind=f64) function, public sll_m_fft::sll_f_fft_get_mode_r2c_1d ( type(sll_t_fft), intent(in)  plan,
real(kind=f64), dimension(0:), intent(in)  data,
integer(kind=i32), intent(in)  k 
)

Function to reconstruct a complex FFT mode from the data of a r2r transform.

Parameters
[in]planFFT planner object
[in]dataReal FFT coefficients (halfcomplex format)
[in]kMode to be extracted (0 <= k <= N-1)
Returns
Complex coefficient of kth mode

Definition at line 311 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_deallocate_aligned_complex()

subroutine, public sll_m_fft::sll_s_fft_deallocate_aligned_complex ( complex(c_double_complex), dimension(:), pointer  data)

Function to deallocate an aligned complex array.

Parameters
dataArray to be deallocated

Definition at line 215 of file sll_m_fft_fftw3.F90.

◆ sll_s_fft_deallocate_aligned_real()

subroutine, public sll_m_fft::sll_s_fft_deallocate_aligned_real ( real(c_double), dimension(:), pointer  data)

Function to deallocate an aligned real array.

Parameters
dataArray to be deallocated

Definition at line 232 of file sll_m_fft_fftw3.F90.

◆ sll_s_fft_exec_c2c_1d()

subroutine, public sll_m_fft::sll_s_fft_exec_c2c_1d ( type(sll_t_fft), intent(in)  plan,
complex(kind=f64), dimension(:), intent(inout)  array_in,
complex(kind=f64), dimension(:), intent(inout)  array_out 
)

Compute fast Fourier transform in complex to complex mode.

Parameters
[in]planFFT planner object
[in,out]array_inComplex data to be Fourier transformed
[in,out]array_outFourier coefficients on output

Definition at line 413 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_exec_c2c_2d()

subroutine, public sll_m_fft::sll_s_fft_exec_c2c_2d ( type(sll_t_fft), intent(in)  plan,
dimension(0:, 0:), intent(inout)  array_in,
dimension(0:, 0:), intent(inout)  array_out 
)

Compute fast Fourier transform in complex to complex mode.

Parameters
[in]planFFT planner object

Definition at line 491 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_exec_c2c_3d()

subroutine, public sll_m_fft::sll_s_fft_exec_c2c_3d ( type(sll_t_fft), intent(in)  plan,
complex(kind=f64), dimension(0:, 0:, 0:), intent(inout)  array_in,
complex(kind=f64), dimension(0:, 0:, 0:), intent(inout)  array_out 
)

Compute fast Fourier transform in complex to complex mode.

Parameters
[in]planFFT planner object
[in,out]array_inComplex data to be Fourier transformed
[in,out]array_outFourier coefficients on output

Definition at line 569 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_exec_c2r_1d()

subroutine, public sll_m_fft::sll_s_fft_exec_c2r_1d ( type(sll_t_fft), intent(in)  plan,
complex(kind=f64), dimension(:), intent(inout)  array_in,
real(kind=f64), dimension(:), intent(inout)  array_out 
)

Compute fast Fourier transform in complex to real mode.

Parameters
[in]planFFT planner objece
[in,out]array_inComplex Fourier coefficient to be transformed back
[in,out]array_outReal result of Fourier transform

Definition at line 914 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_exec_c2r_2d()

subroutine, public sll_m_fft::sll_s_fft_exec_c2r_2d ( type(sll_t_fft), intent(in)  plan,
complex(kind=f64), dimension(1:, 1:), intent(inout)  array_in,
real(kind=f64), dimension(1:, 1:), intent(out)  array_out 
)

Compute fast Fourier transform in complex to real mode.

Parameters
[in]planFFT planner object
[in,out]array_inComplex Fourier coefficient to be transformed back
[out]array_outReal output of Fourier transform

Definition at line 982 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_exec_r2c_1d()

subroutine, public sll_m_fft::sll_s_fft_exec_r2c_1d ( type(sll_t_fft), intent(in)  plan,
real(kind=f64), dimension(:), intent(inout)  array_in,
complex(kind=f64), dimension(:), intent(out)  array_out 
)

Compute fast Fourier transform in real to complex mode.

Parameters
[in]planFFT planner object
[in,out]array_inReal input data to be Fourier transformed
[out]array_outComplex Fourier mode (only first half due to symmetry)

Definition at line 754 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_exec_r2c_2d()

subroutine, public sll_m_fft::sll_s_fft_exec_r2c_2d ( type(sll_t_fft), intent(in)  plan,
real(kind=f64), dimension(:, :), intent(inout)  array_in,
complex(kind=f64), dimension(:, :), intent(out)  array_out 
)

Compute fast Fourier transform in real to complex mode.

Parameters
[in]planFFT planner object
[in,out]array_inReal input data to be Fourier transformed
[out]array_outComplex Fourier coefficients (only half part along first dimension due to symmetry)

Definition at line 824 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_exec_r2r_1d()

subroutine, public sll_m_fft::sll_s_fft_exec_r2r_1d ( type(sll_t_fft), intent(inout)  plan,
real(kind=f64), dimension(:), intent(inout)  array_in,
real(kind=f64), dimension(:), intent(inout)  array_out 
)

Compute fast Fourier transform in real to real mode.

Parameters
[in,out]planFFT planner object
[in,out]array_inReal data to be Fourier transformed
[in,out]array_outFourier coefficients in real form (sin/cos coefficients)

Definition at line 660 of file sll_m_fft_fftw3.F90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sll_s_fft_free()

subroutine, public sll_m_fft::sll_s_fft_free ( type(sll_t_fft), intent(inout)  plan)

Delete a plan.

Parameters
[in,out]planFFT planner object

Definition at line 1006 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_get_k_list_c2c_1d()

subroutine, public sll_m_fft::sll_s_fft_get_k_list_c2c_1d ( type(sll_t_fft), intent(in)  plan,
integer(kind=i32), dimension(0:), intent(out)  k_list 
)

Get a list of the k modes in a c2c FFT k_list = [0, 1, 2, ..., n/2, (n/2+1)-n, ..., -2, -1].

Parameters
[in]planFFT planner object
[out]k_listArray with all the k modes

Definition at line 250 of file sll_m_fft_fftw3.F90.

◆ sll_s_fft_get_k_list_r2c_1d()

subroutine, public sll_m_fft::sll_s_fft_get_k_list_r2c_1d ( type(sll_t_fft), intent(in)  plan,
integer(kind=i32), dimension(0:), intent(out)  k_list 
)

Get a list of the k modes in an r2c FFT k_list = [0, 1, 2, ..., n/2].

Parameters
[in]planFFT planner object
[out]k_listArray with all the k modes

Definition at line 292 of file sll_m_fft_fftw3.F90.

◆ sll_s_fft_get_k_list_r2r_1d()

subroutine, public sll_m_fft::sll_s_fft_get_k_list_r2r_1d ( type(sll_t_fft), intent(in)  plan,
integer(kind=i32), dimension(0:), intent(out)  k_list 
)

Get a list of the k modes in an r2r FFT k_list = [0, 1, 2, ..., n/2, (n+1)/2-1, ..., 2, 1].

Parameters
[in]planFFT planner object
[out]k_listArray with all the k modes

Definition at line 271 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_init_c2c_1d()

subroutine, public sll_m_fft::sll_s_fft_init_c2c_1d ( type(sll_t_fft plan,
integer(kind=i32), intent(in)  nx,
complex(kind=f64), dimension(:), intent(inout)  array_in,
complex(kind=f64), dimension(:), intent(inout)  array_out,
integer(kind=i32), intent(in)  direction,
logical, intent(in), optional  normalized,
logical, intent(in), optional  aligned,
integer(kind=i32), intent(in), optional  optimization 
)

Create new 1d complex to complex plan.

Parameters
planFFT planner object
[in]nxNumber of points
[in,out]array_in(Typical) input array (gets overwritten for certain options)
[in,out]array_out(Typical) output array (gets overwritten for certain options)
[in]directionDirection of the FFT (sll_p_fft_forward or sll_p_fft_backward)
[in]normalizedFlag to decide if FFT should be normalized by 1/N (default: FALSE)
[in]alignedFlag to decide if FFT routine can assume data alignment (default: FALSE). Note that you need to call an aligned initialization if you want to set this option to TRUE.
[in]optimizationPlanning-rigor flag for FFTW. Possible values sll_p_fft_estimate, sll_p_fft_measure, sll_p_fft_patient, sll_p_fft_exhaustive, sll_p_fft_wisdom_only. (default: sll_p_fft_estimate).

Definition at line 361 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_init_c2c_2d()

subroutine, public sll_m_fft::sll_s_fft_init_c2c_2d ( type(sll_t_fft), intent(out)  plan,
integer(kind=i32), intent(in)  nx,
integer(kind=i32), intent(in)  ny,
complex(kind=f64), dimension(0:, 0:), intent(inout)  array_in,
complex(kind=f64), dimension(0:, 0:), intent(inout)  array_out,
integer(kind=i32), intent(in)  direction,
logical, intent(in), optional  normalized,
logical, intent(in), optional  aligned,
integer(kind=i32), intent(in), optional  optimization 
)

Create new 2d complex to complex plan.

Parameters
[out]planinitialized planner object
[in]nxNumber of points along first dimension
[in]nyNumber of points along second dimension
[in,out]array_in(Typical) input array (gets overwritten for certain options)
[in,out]array_out(Typical) output array (gets overwritten for certain options)
[in]directionDirection of the FFT (sll_p_fft_forward or sll_p_fft_backward)
[in]normalizedFlag to decide if FFT should be normalized by 1/N (default: FALSE)
[in]alignedFlag to decide if FFT routine can assume data alignment (default: FALSE). Not that you need to call an aligned initialization if you want to set this option to TRUE.
[in]optimizationPlanning-rigor flag for FFTW. Possible values sll_p_fft_estimate, sll_p_fft_measure, sll_p_fft_patient, sll_p_fft_exhaustive, sll_p_fft_wisdom_only. (default: sll_p_fft_estimate). Note that you need to

Definition at line 439 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_init_c2c_3d()

subroutine, public sll_m_fft::sll_s_fft_init_c2c_3d ( type(sll_t_fft), intent(out)  plan,
intent(in)  nx,
intent(in)  ny,
intent(in)  nz,
dimension(0:, 0:, 0:), intent(inout)  array_in,
dimension(0:, 0:, 0:), intent(inout)  array_out,
intent(in)  direction,
logical, intent(in), optional  normalized,
logical, intent(in), optional  aligned,
intent(in)  optimization 
)

Create new 3d complex to complex plan.

Parameters
[out]planinitialized planner object
[in]normalizedFlag to decide if FFT should be normalized by 1/N (default: FALSE)
[in]alignedFlag to decide if FFT routine can assume data alignment (default: FALSE). Not that you need to call an aligned initialization if you want to set this option to TRUE.

Definition at line 516 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_init_c2r_1d()

subroutine, public sll_m_fft::sll_s_fft_init_c2r_1d ( type(sll_t_fft), intent(out)  plan,
integer(kind=i32), intent(in)  nx,
complex(kind=f64), dimension(:), intent(inout)  array_in,
real(kind=f64), dimension(:), intent(out)  array_out,
logical, intent(in), optional  normalized,
logical, intent(in), optional  aligned,
integer(kind=i32), intent(in), optional  optimization 
)

Create new 1d complex to real plan for backward FFT.

Parameters
[out]planFFT planner object
[in]nxNumber of points
[in,out]array_in(Typical) input array (gets overwritten for certain options)
[out]array_out(Typical) output array (gets overwritten for certain options)
[in]normalizedFlag to decide if FFT should be normalized by 1/N (default: FALSE)
[in]alignedFlag to decide if FFT routine can assume data alignment (default: FALSE). Not that you need to call an aligned initialization if you want to set this option to TRUE.
[in]optimizationPlanning-rigor flag for FFTW. Possible values sll_p_fft_estimate, sll_p_fft_measure, sll_p_fft_patient, sll_p_fft_exhaustive, sll_p_fft_wisdom_only. (default: sll_p_fft_estimate). Note that you need to

Definition at line 866 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_init_c2r_2d()

subroutine, public sll_m_fft::sll_s_fft_init_c2r_2d ( type(sll_t_fft), intent(out)  plan,
integer(kind=i32), intent(in)  nx,
integer(kind=i32), intent(in)  ny,
complex(kind=f64), dimension(:, :), intent(inout)  array_in,
real(kind=f64), dimension(:, :), intent(out)  array_out,
logical, intent(in), optional  normalized,
logical, intent(in), optional  aligned,
integer(kind=i32), intent(in), optional  optimization 
)

Create new 2d real to complex plan for backward FFT.

Parameters
[out]planFFT planner object
[in]nxNumber of point along first dimension
[in]nyNumber of points along second dimension
[in,out]array_in(Typical) input array (gets overwritten for certain options)
[out]array_out(Typical) output array (gets overwritten for certain options)
[in]normalizedFlag to decide if FFT should be normalized by 1/N (default: FALSE)
[in]alignedFlag to decide if FFT routine can assume data alignment (default: FALSE). Not that you need to call an aligned initialization if you want to set this option to TRUE.
[in]optimizationPlanning-rigor flag for FFTW. Possible values sll_p_fft_estimate, sll_p_fft_measure, sll_p_fft_patient, sll_p_fft_exhaustive, sll_p_fft_wisdom_only. (default: sll_p_fft_estimate). Note that you need to

Definition at line 934 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_init_r2c_1d()

subroutine, public sll_m_fft::sll_s_fft_init_r2c_1d ( type(sll_t_fft), intent(out)  plan,
integer(kind=i32), intent(in)  nx,
real(kind=f64), dimension(:), intent(inout)  array_in,
complex(kind=f64), dimension(:), intent(out)  array_out,
logical, intent(in), optional  normalized,
logical, intent(in), optional  aligned,
integer(kind=i32), intent(in), optional  optimization 
)

Create new 1d real to complex plan for forward FFT.

Parameters
[out]planFFT planner object
[in]nxNumber of points
[in,out]array_in(Typical) input array (gets overwritten for certain options)
[out]array_out(Typical) output array (gets overwritten for certain options)
[in]normalizedFlag to decide if FFT should be normalized by 1/N (default: FALSE)
[in]alignedFlag to decide if FFT routine can assume data alignment (default: FALSE). Not that you need to call an aligned initialization if you want to set this option to TRUE.
[in]optimizationPlanning-rigor flag for FFTW. Possible values sll_p_fft_estimate, sll_p_fft_measure, sll_p_fft_patient, sll_p_fft_exhaustive, sll_p_fft_wisdom_only. (default: sll_p_fft_estimate). Note that you need to

Definition at line 709 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_init_r2c_2d()

subroutine, public sll_m_fft::sll_s_fft_init_r2c_2d ( type(sll_t_fft), intent(out)  plan,
integer(kind=i32), intent(in)  nx,
integer(kind=i32), intent(in)  ny,
real(kind=f64), dimension(:, :), intent(inout)  array_in,
complex(kind=f64), dimension(:, :), intent(out)  array_out,
logical, intent(in), optional  normalized,
logical, intent(in), optional  aligned,
integer(kind=i32), intent(in), optional  optimization 
)

Create new 2d complex to real plan for forward FFT.

Parameters
[out]planFFT planner object
[in]nxNumber of points along first dimension
[in]nyNumber of points along second dimension
[in,out]array_in(Typical) input array (gets overwritten for certain options)
[out]array_out(Typical) output array (gets overwritten for certain options)
[in]normalizedFlag to decide if FFT should be normalized by 1/N (default: FALSE)
[in]alignedFlag to decide if FFT routine can assume data alignment (default: FALSE). Not that you need to call an aligned initialization if you want to set this option to TRUE.
[in]optimizationPlanning-rigor flag for FFTW. Possible values sll_p_fft_estimate, sll_p_fft_measure, sll_p_fft_patient, sll_p_fft_exhaustive, sll_p_fft_wisdom_only. (default: sll_p_fft_estimate). Note that you need to

Definition at line 774 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_init_r2r_1d()

subroutine, public sll_m_fft::sll_s_fft_init_r2r_1d ( type(sll_t_fft), intent(out)  plan,
integer(kind=i32), intent(in)  nx,
real(kind=f64), dimension(:), intent(inout)  array_in,
real(kind=f64), dimension(:), intent(inout)  array_out,
integer(kind=i32), intent(in)  direction,
logical, intent(in), optional  normalized,
logical, intent(in), optional  aligned,
integer(kind=i32), intent(in), optional  optimization 
)

Create new 1d real to real plan.

Parameters
[out]planFFT planner object
[in]nxNumber of points
[in,out]array_in(Typical) input array (gets overwritten for certain options)
[in,out]array_out(Typical) output array (gets overwritten for certain options)
[in]directionDirection of the FFT (sll_p_fft_forward or sll_p_fft_backward)
[in]normalizedFlag to decide if FFT should be normalized by 1/N (default: FALSE)
[in]alignedFlag to decide if FFT routine can assume data alignment (default: FALSE). Not that you need to call an aligned initialization if you want to set this option to TRUE.
[in]optimizationPlanning-rigor flag for FFTW. Possible values sll_p_fft_estimate, sll_p_fft_measure, sll_p_fft_patient, sll_p_fft_exhaustive, sll_p_fft_wisdom_only. (default: sll_p_fft_estimate). Note that you need to

Definition at line 594 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

◆ sll_s_fft_print_defaultlib()

subroutine, public sll_m_fft::sll_s_fft_print_defaultlib

Function to print the FFT library behind the interface.

Definition at line 170 of file sll_m_fft_fftw3.F90.

◆ sll_s_fft_set_mode_c2r_1d()

subroutine, public sll_m_fft::sll_s_fft_set_mode_c2r_1d ( type(sll_t_fft), intent(in)  plan,
real(kind=f64), dimension(0:), intent(inout)  data,
complex(kind=f64), intent(in)  ck,
integer(kind=i32), intent(in)  k 
)

Subroutine to set a complex mode to the real representation of r2r.

Parameters
[in]planFFT planner object
[in,out]dataReal FFT coefficients (halfcomplex format)
[in]ckComplex coefficient of kth mode
[in]kMode to be set (0 <= k <= N-1)

Definition at line 333 of file sll_m_fft_fftw3.F90.

Here is the caller graph for this function:

Variable Documentation

◆ p_fftw_c2c_1d

integer, parameter p_fftw_c2c_1d = 0
private

Definition at line 153 of file sll_m_fft_fftw3.F90.

◆ p_fftw_c2c_2d

integer, parameter p_fftw_c2c_2d = 4
private

Definition at line 157 of file sll_m_fft_fftw3.F90.

◆ p_fftw_c2c_3d

integer, parameter p_fftw_c2c_3d = 8
private

Definition at line 161 of file sll_m_fft_fftw3.F90.

◆ p_fftw_c2r_1d

integer, parameter p_fftw_c2r_1d = 3
private

Definition at line 156 of file sll_m_fft_fftw3.F90.

◆ p_fftw_c2r_2d

integer, parameter p_fftw_c2r_2d = 7
private

Definition at line 160 of file sll_m_fft_fftw3.F90.

◆ p_fftw_r2c_1d

integer, parameter p_fftw_r2c_1d = 2
private

Definition at line 155 of file sll_m_fft_fftw3.F90.

◆ p_fftw_r2c_2d

integer, parameter p_fftw_r2c_2d = 6
private

Definition at line 159 of file sll_m_fft_fftw3.F90.

◆ p_fftw_r2r_1d

integer, parameter p_fftw_r2r_1d = 1
private

Definition at line 154 of file sll_m_fft_fftw3.F90.

◆ p_fftw_r2r_2d

integer, parameter p_fftw_r2r_2d = 5
private

Definition at line 158 of file sll_m_fft_fftw3.F90.

◆ sll_p_fft_backward

integer, parameter, public sll_p_fft_backward = FFTW_BACKWARD

Flag to specify backward FFT (positive sign) [value 1].

Definition at line 138 of file sll_m_fft_fftw3.F90.

◆ sll_p_fft_estimate

integer(c_int), parameter, public sll_p_fft_estimate = FFTW_ESTIMATE

FFTW planning-rigor flag FFTW_ESTIMATE (simple heuristic for planer) [value 64]. This is our default value.

Definition at line 143 of file sll_m_fft_fftw3.F90.

◆ sll_p_fft_exhaustive

integer(c_int), parameter, public sll_p_fft_exhaustive = FFTW_EXHAUSTIVE

FFTW planning-rigor flag FFTW_EXHAUSTIVE (more optimization than PATIENT) NOTE: planner overwrites the input array during planning [value 8].

Definition at line 144 of file sll_m_fft_fftw3.F90.

◆ sll_p_fft_forward

integer, parameter, public sll_p_fft_forward = FFTW_FORWARD

Flag to specify forward FFT (negative sign) [value -1].

Definition at line 137 of file sll_m_fft_fftw3.F90.

◆ sll_p_fft_measure

integer(c_int), parameter, public sll_p_fft_measure = FFTW_MEASURE

FFTW planning-rigor flag FFTW_MEASURE (optimized plan) NOTE: planner overwrites the input array during planning [value 0].

Definition at line 141 of file sll_m_fft_fftw3.F90.

◆ sll_p_fft_patient

integer(c_int), parameter, public sll_p_fft_patient = FFTW_PATIENT

FFTW planning-rigor flag FFTW_PATIENT (more optimizaton than MEASURE) NOTE: planner overwrites the input array during planning [value 32].

Definition at line 142 of file sll_m_fft_fftw3.F90.

◆ sll_p_fft_wisdom_only

integer(c_int), parameter, public sll_p_fft_wisdom_only = FFTW_WISDOM_ONLY

Definition at line 145 of file sll_m_fft_fftw3.F90.

    Report Typos and Errors