|
| | R2C (Box const &domain, Info const &info=Info{}) |
| | Constructor.
|
| |
| | R2C (std::array< int, 3 > const &domain_size, Info const &info=Info{}) |
| | Constructor.
|
| |
| | ~R2C () |
| |
| | R2C (R2C const &)=delete |
| |
| | R2C (R2C &&)=delete |
| |
| R2C & | operator= (R2C const &)=delete |
| |
| R2C & | operator= (R2C &&)=delete |
| |
| void | setLocalDomain (std::array< int, 3 > const &local_start, std::array< int, 3 > const &local_size) |
| | Set local domain.
|
| |
| std::pair< std::array< int, 3 >, std::array< int, 3 > > | getLocalDomain () const |
| | Get local domain.
|
| |
| void | setLocalSpectralDomain (std::array< int, 3 > const &local_start, std::array< int, 3 > const &local_size) |
| | Set local spectral domain.
|
| |
| std::pair< std::array< int, 3 >, std::array< int, 3 > > | getLocalSpectralDomain () const |
| | Get local spectral domain.
|
| |
| 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, int incomp=0, int outcomp=0) |
| | Forward and then backward transform.
|
| |
| template<Direction DIR = D, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > = 0> |
| void | forward (MF const &inmf, int incomp=0) |
| | Forward transform.
|
| |
| template<Direction DIR = D, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > = 0> |
| void | forward (MF const &inmf, cMF &outmf, int incomp=0, int outcomp=0) |
| | Forward transform.
|
| |
| template<typename RT , typename CT , Direction DIR = D, bool CP = C, std::enable_if_t<(DIR==Direction::forward||DIR==Direction::both) &&((sizeof(RT) *2==sizeof(CT) &&!CP)||(sizeof(RT)==sizeof(CT) &&CP)), int > = 0> |
| void | forward (RT const *in, CT *out) |
| | Forward transform.
|
| |
| template<Direction DIR = D, std::enable_if_t< DIR==Direction::both, int > = 0> |
| void | backward (MF &outmf, int outcomp=0) |
| | Backward transform.
|
| |
| template<Direction DIR = D, std::enable_if_t< DIR==Direction::backward||DIR==Direction::both, int > = 0> |
| void | backward (cMF const &inmf, MF &outmf, int incomp=0, int outcomp=0) |
| | Backward transform.
|
| |
| template<typename CT , typename RT , Direction DIR = D, bool CP = C, std::enable_if_t<(DIR==Direction::backward||DIR==Direction::both) &&((sizeof(RT) *2==sizeof(CT) &&!CP)||(sizeof(RT)==sizeof(CT) &&CP)), int > = 0> |
| void | backward (CT const *in, RT *out) |
| | Backward transform.
|
| |
| T | scalingFactor () const |
| |
| template<Direction DIR = D, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > = 0> |
| std::pair< cMF *, IntVect > | getSpectralData () const |
| | Get the internal spectral data.
|
| |
| std::pair< BoxArray, DistributionMapping > | getSpectralDataLayout () const |
| | Get BoxArray and DistributionMapping for spectral data.
|
| |
| template<typename F > |
| void | post_forward_doit_0 (F const &post_forward) |
| |
| template<typename F > |
| void | post_forward_doit_1 (F const &post_forward) |
| |
| template<Direction DIR, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > > |
| std::pair< typename R2C< T, D, C >::cMF *, IntVect > | getSpectralData () const |
| |
template<typename T = Real,
FFT::Direction D = FFT::Direction::both, bool C = false>
class amrex::FFT::R2C< T, D, C >
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.
The arrays are assumed to be in column-major, which is different from FFTW's row-major layout. Because the complex domain data have the Hermitian symmetry, only half of the data in the complex domain are stored. If the real domain size is nx * ny * nz, the complex domain's size will be (nx/2+1) * ny * nz.
For more details, we refer the users to https://amrex-codes.github.io/amrex/docs_html/FFT_Chapter.html.
| std::pair< std::array< int, 3 >, std::array< int, 3 > > amrex::FFT::R2C< T, D, C >::getLocalSpectralDomain |
( |
| ) |
const |
Get local spectral domain.
This function returns the domain decomposition chosen by AMReX for the complex data spectral domain. The returned pair contains the local starting indices and the local domain size.
This is needed, only if one uses the raw pointer interfaces, not the amrex::MulitFab interfaces. Only one of the functions, setLocalSpectralDomain and getLocalSpectralDomain, should be called. Note that one could use this function together with setLocalDomain. That is the user is allowed to choose their own real domain decomposition, while let AMReX choose the spectral data domain decomposition. Also note that the entire spectral domain has the size of (nx+1)/2 * ny * nz, if the real domain is nx * ny * nz.
| void amrex::FFT::R2C< T, D, C >::setLocalDomain |
( |
std::array< int, 3 > const & |
local_start, |
|
|
std::array< int, 3 > const & |
local_size |
|
) |
| |
Set local domain.
This is needed, only if one uses the raw pointer interfaces, not the amrex::MulitFab interfaces. This may contain collective MPI calls. So all processes even if their local size is zero should call. This function informs AMReX the domain decomposition chosen by the user for the forward domain (i.e., real data domain). There is no constraint on the domain decomposition strategy. One can do 1D, 2D or 3D domain decomposition. Alternatively, one could also let AMReX choose for you by calling getLocalDomain. The latter could potentially reduce data communication.
Again, this is needed, only if one uses the raw pointer interfaces. Only one of the functions, setLocalDomain and getLocalDomain, should be called.
This should only be called once unless the domain decomposition changes.
local_start starting indices of the local domain local_size size of the local domain
| void amrex::FFT::R2C< T, D, C >::setLocalSpectralDomain |
( |
std::array< int, 3 > const & |
local_start, |
|
|
std::array< int, 3 > const & |
local_size |
|
) |
| |
Set local spectral domain.
This is needed, only if one uses the raw pointer interfaces, not the amrex::MulitFab interfaces. This may contain collective MPI calls. So all processes even if their local size is zero should call. This function informs AMReX the domain decomposition chosen by the user for the complex data domain. There is no constraint on the domain decomposition strategy. One can do 1D, 2D or 3D domain decomposition. Alternatively, one could also let AMReX choose for you by calling getLocalSpectralDomain. The latter could potentially reduce data communication.
Again, this is needed, only if one uses the raw pointer interfaces. Only one of the functions, setLocalSpectralDomain and getLocalSpectralDomain, should be called. Note that one could use this function together with getLocalDomain. That is the user is allowed to choose their own spectral domain decomposition, while let AMReX choose the real data domain decomposition. Also note that the entire spectral domain has the size of (nx+1)/2 * ny * nz, if the real domain is nx * ny * nz.
This should only be called once unless the domain decomposition changes.
local_start starting indices of the local domain local_size size of the local domain