Local Discrete Fourier Transform. More...
#include <AMReX_FFT_LocalR2C.H>
Public Member Functions | |
LocalR2C (IntVectND< M > const &fft_size, T *p_fwd=nullptr, GpuComplex< T > *p_bwd=nullptr, bool cache_plan=true) | |
Constructor. | |
~LocalR2C () | |
LocalR2C ()=default | |
LocalR2C (LocalR2C &&) noexcept | |
LocalR2C & | operator= (LocalR2C &&) noexcept |
LocalR2C (LocalR2C const &)=delete | |
LocalR2C & | operator= (LocalR2C const &)=delete |
template<Direction DIR = D, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > = 0> | |
void | forward (T const *indata, GpuComplex< T > *outdata) |
Forward transform. | |
void | clear () |
template<Direction DIR = D, std::enable_if_t< DIR==Direction::backward||DIR==Direction::both, int > = 0> | |
void | backward (GpuComplex< T > const *indata, T *outdata) |
Backward transform. | |
T | scalingFactor () const |
IntVectND< M > const & | spectralSize () const |
Spectral domain size. | |
Private Attributes | |
Plan< T > | m_fft_fwd |
Plan< T > | m_fft_bwd |
T * | m_p_fwd = nullptr |
GpuComplex< T > * | m_p_bwd = nullptr |
IntVectND< M > | m_real_size |
IntVectND< M > | m_spectral_size |
bool | m_cache_plan = false |
Local Discrete Fourier Transform.
This class supports Fourier transforms between real and complex data. The name R2C indicates that the forward transform converts real data to complex data, while the backward transform converts complex data to real data. It should be noted that both directions of transformation are supported, not just from real to complex. The scaling follows the FFTW convention, where applying the forward transform followed by the backward transform scales the original data by the size of the input array.
For more details, we refer the users to https://amrex-codes.github.io/amrex/docs_html/FFT_Chapter.html.
|
explicit |
Constructor.
Given the diverse interfaces of FFT libraries we use, this constructo has a number of optional arguments.
The user can provide the data pointers to the constructor. They are only needed by FFTW because its plan creation requires the input and output arrays. If they are null, we will delay the plan creation for FFTW until the forward or backward function is called.
The cache_plan option is only used when we use cufft, rocfft and onemkl, but not FFTW.
fft_size | The forward domain size (i.e., the domain of the real data) |
p_fwd | Forward domain data pointer (optional) |
p_bwd | Backward domain data pointer (optional) |
cache_plan | Try to cache the plan or not (optionl) |
amrex::FFT::LocalR2C< T, D, M >::~LocalR2C | ( | ) |
|
default |
|
noexcept |
|
delete |
void amrex::FFT::LocalR2C< T, D, M >::backward | ( | GpuComplex< T > const * | indata, |
T * | outdata | ||
) |
Backward transform.
This function is not available when this class template is instantiated for forward-only transform. For GPUs, this function is synchronous on the host.
indata | input data |
outdata | output data |
void amrex::FFT::LocalR2C< T, D, M >::clear | ( | ) |
void amrex::FFT::LocalR2C< T, D, M >::forward | ( | T const * | indata, |
GpuComplex< T > * | outdata | ||
) |
Forward transform.
This function is not available when this class template is instantiated for backward-only transform. For GPUs, this function is synchronous on the host.
indata | input data |
outdata | output data |
|
noexcept |
|
delete |
T amrex::FFT::LocalR2C< T, D, M >::scalingFactor | ( | ) | const |
Scaling factor. If the data goes through forward and then backward, the result multiplied by the scaling factor is equal to the original data.
|
inline |
Spectral domain size.
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |