|
| R2C (Box const &domain, Info const &info=Info{}) |
| Constructor. More...
|
|
| ~R2C () |
|
| R2C (R2C const &)=delete |
|
| R2C (R2C &&)=delete |
|
R2C & | operator= (R2C const &)=delete |
|
R2C & | operator= (R2C &&)=delete |
|
template<typename F , Direction DIR = D, std::enable_if_t< DIR==Direction::both, int > = 0> |
void | forwardThenBackward (MF const &inmf, MF &outmf, F const &post_forward) |
| Forward and then backward transform. More...
|
|
template<Direction DIR = D, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > = 0> |
void | forward (MF const &inmf) |
| Forward transform. More...
|
|
template<Direction DIR = D, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > = 0> |
void | forward (MF const &inmf, cMF &outmf) |
| Forward transform. More...
|
|
template<Direction DIR = D, std::enable_if_t< DIR==Direction::both, int > = 0> |
void | backward (MF &outmf) |
| Backward transform. More...
|
|
template<Direction DIR = D, std::enable_if_t< DIR==Direction::backward||DIR==Direction::both, int > = 0> |
void | backward (cMF const &inmf, MF &outmf) |
| Backward transform. More...
|
|
T | scalingFactor () const |
|
template<Direction DIR = D, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > = 0> |
std::pair< cMF *, IntVect > | getSpectralData () |
| Get the internal spectral data. More...
|
|
std::pair< BoxArray, DistributionMapping > | getSpectralDataLayout () const |
| Get BoxArray and DistributionMapping for spectral data. More...
|
|
template<int Depth, typename F > |
void | post_forward_doit (F const &post_forward) |
|
template<Direction DIR, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > > |
std::pair< typename R2C< T, D, S >::cMF *, IntVect > | getSpectralData () |
|
template<typename T = Real, FFT::Direction D = FFT::Direction::both, FFT::DomainStrategy S = FFT::DomainStrategy::slab>
class amrex::FFT::R2C< T, D, S >
Parallel 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.
template<typename T , Direction D, DomainStrategy S>
template<Direction DIR, std::enable_if_t< DIR==Direction::backward||DIR==Direction::both, int > >
Backward transform.
This function is not available when this class template is instantiated for forward-only transform.
- Parameters
-
template<typename T , Direction D, DomainStrategy S>
template<Direction DIR, std::enable_if_t< DIR==Direction::both, int > >
Backward transform.
This function is available only when this class template is instantiated for transforms in both directions.
- Parameters
-
template<typename T , Direction D, DomainStrategy S>
template<Direction DIR, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > >
Forward transform.
The output is stored in this object's internal data. This function is not available when this class template is instantiated for backward-only transform.
- Parameters
-
template<typename T = Real, FFT::Direction D = FFT::Direction::both, FFT::DomainStrategy S = FFT::DomainStrategy::slab>
template<typename F , Direction DIR = D, std::enable_if_t< DIR==Direction::both, int > = 0>
void amrex::FFT::R2C< T, D, S >::forwardThenBackward |
( |
MF const & |
inmf, |
|
|
MF & |
outmf, |
|
|
F const & |
post_forward |
|
) |
| |
|
inline |
Forward and then backward transform.
This function is available only when this class template is instantiated for transforms in both directions. It's more efficient than calling the forward function that stores the spectral data in a caller provided container followed by the backward function, because this can avoid parallel communication between the internal data and the caller's data container.
- Parameters
-
inmf | input data in MultiFab or FabArray<BaseFab<float>> |
outmf | output data in MultiFab or FabArray<BaseFab<float>> |
post_forward | a callable object for processing the post-forward data before the backward transform. Its interface is (int,int,int,GpuComplex<T>&) , where the integers are indices in the spectral space, and the reference to the complex number allows for the modification of the spectral data at that location. |
template<typename T = Real, FFT::Direction D = FFT::Direction::both, FFT::DomainStrategy S = FFT::DomainStrategy::slab>
template<Direction DIR = D, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > = 0>
Get the internal spectral data.
This function is not available when this class template is instantiated for backward-only transform. For performance reasons, the returned data array does not have the usual ordering of (x,y,z)
. The order is specified in the second part of the return value.