SF_FFT
Description
SciFortran module for Fourier transform
Quick access
- Routines:
cfft_1d_backward(),cfft_1d_ex(),cfft_1d_forward(),cfft_1d_ishift(),cfft_1d_shift(),cfft_2d_backward(),cfft_2d_forward(),cfft_nd_backward(),cfft_nd_forward(),cosft(),cosftn(),cost_1d_backward(),cost_1d_forward(),cost_nd_backward(),cost_nd_forward(),fft(),fft2(),fft_farray(),fft_fmax(),fft_signal(),fft_tarray(),fft_tmax(),fftex(),fftn(),fftshift(),ft_direct(),ft_inverse(),icosft(),icosftn(),ifft(),ifft2(),ifft_signal(),ifftn(),ifftshift(),isinft(),isinftn(),itfft(),rfft_1d_backward(),rfft_1d_ex(),rfft_1d_forward(),rfft_1d_ishift(),rfft_1d_shift(),rfft_2d_backward(),rfft_2d_forward(),rfft_nd_backward(),rfft_nd_forward(),sinft(),sinftn(),sint_1d_backward(),sint_1d_forward(),sint_nd_backward(),sint_nd_forward(),tfft()
Used modules
sf_integrate: SciFortran module for function integrationsf_arrays: SciFortran module for array creation and manipulationsf_constants: SciFortran module for physical and mathematical constants
Subroutines and functions
- interface sf_fft_fftpack/ft_direct(ft, t, w)
This function evaluates the direct Fourier transform of a discretized function from time to frequency domain. Takes as input a discretized function \(ft_{i}\) on a set of time points \(t_{i}\) and a set of frequencies \(w_{j}\) and returns a discretized function \(fw_{j}\) =
simps\((ft_{i} \cdot e^{-i 2\pi w_{j}t_{i}}, t_{1}, t_{N})\) , whereN=size(t)- Parameters:
ft (•) [real/complex, in,required] – Dicretized function to transform (time domain)
t (size(ft)) [real, in] – Discretized time points
w (•) [real, in] – Discretized frequency points
- Result:
fw (size(w)) [real/complex] – Discretized Fourier-transform function (frequency domain)
- interface sf_fft_fftpack/ft_inverse(fw, t, w)
This function evaluates the inverse Fourier transform of a discretized function from frequency to time domain. Takes as input a discretized function \(fw_{i}\) on a set of time points \(w_{i}\) and a set of frequencies \(t_{j}\) and returns a discretized function \(ft_{j}\) =
simps\((fw_{i} \cdot e^{i 2\pi t_{j}w_{i}}, w_{1}, w_{N})\) , whereN=size(w)- Parameters:
fw (•) [real/complex, in,required] – Dicretized function to transform (frequency domain)
t (•) [real, in] – Discretized time points
w (size(fw)) [real, in] – Discretized frequency points
- Result:
ft (size(t)) [real/complex] – Discretized Fourier-transform function (time domain)
- interface sf_fft_fftpack/fft_signal(ft, dt)
This function evaluates the fast Fourier transform of a discretized function
ft. It returnsfw=dt·tfft(ft)- Parameters:
ft (•) [real/complex] – Time-domain function
dt [real] – Time step
- Result:
fw (size(ft)) [real/complex] – Frequency-domain function
- interface sf_fft_fftpack/ifft_signal(fw, dt)
This function evaluates the inverse fast Fourier transform of a discretized function
fw. It returnsft=itfft(fw) /ft- Parameters:
fw (•) [real/complex] – Frequency-domain function
dt [real] – Time step
- Result:
ft (size(fw)) [real/complex] – Time-domain function
- interface sf_fft_fftpack/tfft(func_in[, func_out])
- Parameters:
func_in (•) [real/complex, in,required]
- Options:
func_out (size(func_in)) [real/complex]
- interface sf_fft_fftpack/itfft(func_in[, func_out])
- Parameters:
func_in (•) [real/complex, in,required]
- Options:
func_out (size(func_in)) [real/complex]
- interface sf_fft_fftpack/fft(func)
This subroutine evaluates the forward 1-dimensional FFT of an array of lenght
Nusing the FFTPACK 5.1 routinesrfft1i+rfft1ffor the real case andcfft1i+cfft1ffor the complex case. These are called with the parameterslensav=N + int( log(dble(N))/log(2.d0) ) + 4lenwrk=Nlenr=Ninc=1
The subroutine modifies the input array to contain
[y(0), y(1), ..., y(N/2), y(-N/2+1), ..., y(-1)]if N is even[y(0), y(1), ..., y((N-1)/2), y(-(N-1)/2), ..., y(-1)]if N is odd
where
\(y(j) = \sum_{k=0}^{N-1} x(k) \cdot e^{-i 2\pi/N \cdot j \cdot k}\)
for \(j \in [0,N-1]\)
- Parameters:
func (•) [real/complex, inout] – Input/output array of size
N
- interface sf_fft_fftpack/ifft(func)
This subroutine evaluates the backward 1-dimensional FFT of an array of lenght
Nusing the FFTPACK 5.1 routinesrfft1i+rfft1bfor the real case andcfft1i+cfft1bfor the complex case. These are called with the parameterslensav=N + int( log(dble(N))/log(2.d0) ) + 4lenwrk=Nlenr=Ninc=1
The subroutine modifies the input array to contain
[y(0), y(1), ..., y(N/2), y(-N/2+1), ..., y(-1)]if N is even[y(0), y(1), ..., y((N-1)/2), y(-(N-1)/2), ..., y(-1)]if N is odd
where
\(y(j) = \sum_{k=-N/2}^{N/2-1} x(k) \cdot a^{i 2\pi/N \cdot j \cdot k}\)
for \(j \in [0,N-1]\)
- Parameters:
func (•) [real/complex, inout]
- interface sf_fft_fftpack/fft2(func)
- Parameters:
func (•, •) [real/complex, inout]
- interface sf_fft_fftpack/ifft2(func)
- Parameters:
func (•, •) [real/complex, inout]
- interface sf_fft_fftpack/fftn(func, n, lot)
- Parameters:
func (•) [real/complex, inout]
n [integer, in]
lot [integer, in]
- interface sf_fft_fftpack/ifftn(func, n, lot)
- Parameters:
func (•) [real/complex, inout]
n [integer, in]
lot [integer, in]
- interface sf_fft_fftpack/cosft(func)
- Parameters:
func (•) [real, inout]
- interface sf_fft_fftpack/icosft(func)
- Parameters:
func (•) [real, inout]
- interface sf_fft_fftpack/cosftn(func, n, lot)
- Parameters:
func (•) [real, inout]
n [integer, in]
lot [integer, in]
- interface sf_fft_fftpack/icosftn(func, n, lot)
- Parameters:
func (•) [real, inout]
n [integer, in]
lot [integer, in]
- interface sf_fft_fftpack/sinft(func)
- Parameters:
func (•) [real, inout]
- interface sf_fft_fftpack/isinft(func)
- Parameters:
func (•) [real, inout]
- interface sf_fft_fftpack/sinftn(func, n, lot)
- Parameters:
func (•) [real, inout]
n [integer, in]
lot [integer, in]
- interface sf_fft_fftpack/isinftn(func, n, lot)
- Parameters:
func (•) [real, inout]
n [integer, in]
lot [integer, in]
- interface sf_fft_fftpack/fftshift(fin)
- Parameters:
fin (•) [real/complex]
- Result:
fout (size(fin)) [real/complex]
- interface sf_fft_fftpack/ifftshift(fin)
- Parameters:
fin (•) [real/complex]
- Result:
fout (size(fin)) [real/complex]
- interface sf_fft_fftpack/fftex(func)
- Parameters:
func (•) [real/complex]
- subroutine sf_fft_fftpack/rfft_1d_forward(func)
- Parameters:
func (•) [real, inout] – Input/output array of size
N
- subroutine sf_fft_fftpack/cfft_1d_forward(func)
- Parameters:
func (•) [complex, inout]
- subroutine sf_fft_fftpack/rfft_2d_forward(func)
- Parameters:
func (•, •) [real, inout]
- subroutine sf_fft_fftpack/cfft_2d_forward(func)
- Parameters:
func (•, •) [complex, inout]
- subroutine sf_fft_fftpack/rfft_nd_forward(func, n, lot)
- Parameters:
func (•) [real, inout]
n [integer, in]
lot [integer, in]
- subroutine sf_fft_fftpack/cfft_nd_forward(func, n, lot)
- Parameters:
func (•) [complex, inout]
n [integer, in]
lot [integer, in]
- subroutine sf_fft_fftpack/cost_1d_forward(func)
- Parameters:
func (•) [real, inout]
- subroutine sf_fft_fftpack/sint_1d_forward(func)
- Parameters:
func (•) [real, inout]
- subroutine sf_fft_fftpack/cost_nd_forward(func, n, lot)
- Parameters:
func (•) [real, inout]
n [integer, in]
lot [integer, in]
- subroutine sf_fft_fftpack/sint_nd_forward(func, n, lot)
- Parameters:
func (•) [real, inout]
n [integer, in]
lot [integer, in]
- subroutine sf_fft_fftpack/rfft_1d_backward(func)
- Parameters:
func (•) [real, inout]
- subroutine sf_fft_fftpack/cfft_1d_backward(func)
- Parameters:
func (•) [complex, inout]
- subroutine sf_fft_fftpack/rfft_2d_backward(func)
- Parameters:
func (•, •) [real, inout]
- subroutine sf_fft_fftpack/cfft_2d_backward(func)
- Parameters:
func (•, •) [complex, inout]
- subroutine sf_fft_fftpack/rfft_nd_backward(func, n, lot)
- Parameters:
func (•) [real, inout]
n [integer, in]
lot [integer, in]
- subroutine sf_fft_fftpack/cfft_nd_backward(func, n, lot)
- Parameters:
func (•) [complex, inout]
n [integer, in]
lot [integer, in]
- subroutine sf_fft_fftpack/cost_1d_backward(func)
- Parameters:
func (•) [real, inout]
- subroutine sf_fft_fftpack/sint_1d_backward(func)
- Parameters:
func (•) [real, inout]
- subroutine sf_fft_fftpack/cost_nd_backward(func, n, lot)
- Parameters:
func (•) [real, inout]
n [integer, in]
lot [integer, in]
- subroutine sf_fft_fftpack/sint_nd_backward(func, n, lot)
- Parameters:
func (•) [real, inout]
n [integer, in]
lot [integer, in]
- function sf_fft_fftpack/rfft_1d_shift(fin)
- Parameters:
fin (•) [real]
- Result:
fout (size(fin)) [real]
- function sf_fft_fftpack/cfft_1d_shift(fin)
- Parameters:
fin (•) [complex]
- Result:
fout (size(fin)) [complex]
- function sf_fft_fftpack/rfft_1d_ishift(fin)
- Parameters:
fin (•) [real]
- Result:
fout (size(fin)) [real]
- function sf_fft_fftpack/cfft_1d_ishift(fin)
- Parameters:
fin (•) [complex]
- Result:
fout (size(fin)) [complex]
- subroutine sf_fft_fftpack/rfft_1d_ex(func)
- Parameters:
func (•) [real]
- subroutine sf_fft_fftpack/cfft_1d_ex(func)
- Parameters:
func (•) [complex]
- function sf_fft_fftpack/fft_tmax(l, dt)
- Parameters:
l [integer]
dt [real]
- Result:
fft_tmax [real]
- function sf_fft_fftpack/fft_fmax(l, dt)
- Parameters:
l [integer]
dt [real]
- Result:
fft_fmax [real]
- function sf_fft_fftpack/fft_tarray(l, dt)
- Parameters:
l [integer]
dt [real]
- Result:
time (l) [real]
- function sf_fft_fftpack/fft_farray(l, dt[, df])
- Parameters:
l [integer]
dt [real]
- Options:
df [real]
- Result:
freq (l) [real]