Block-Structured AMR Software Framework
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Loading...
Searching...
No Matches
amrex::FFT::R2C< T, D, C > Class Template Reference

Parallel Discrete Fourier Transform. More...

#include <AMReX_FFT_R2C.H>

Public Types

using cMF = FabArray< BaseFab< GpuComplex< T > > >
 
using MF = std::conditional_t< C, cMF, std::conditional_t< std::is_same_v< T, Real >, MultiFab, FabArray< BaseFab< T > > > >
 

Public Member Functions

 R2C (Box const &domain, Info const &info=Info{})
 Constructor.
 
 R2C (std::array< int, AMREX_SPACEDIM > const &domain_size, Info const &info=Info{})
 Constructor.
 
 ~R2C ()
 
 R2C (R2C const &)=delete
 
 R2C (R2C &&)=delete
 
R2Coperator= (R2C const &)=delete
 
R2Coperator= (R2C &&)=delete
 
void setLocalDomain (std::array< int, AMREX_SPACEDIM > const &local_start, std::array< int, AMREX_SPACEDIM > const &local_size)
 Set local domain.
 
std::pair< std::array< int, AMREX_SPACEDIM >, std::array< int, AMREX_SPACEDIM > > getLocalDomain () const
 Get local domain.
 
void setLocalSpectralDomain (std::array< int, AMREX_SPACEDIM > const &local_start, std::array< int, AMREX_SPACEDIM > const &local_size)
 Set local spectral domain.
 
std::pair< std::array< int, AMREX_SPACEDIM >, std::array< int, AMREX_SPACEDIM > > 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.
 
scalingFactor () const
 
template<Direction DIR = D, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > = 0>
std::pair< cMF *, IntVectgetSpectralData () const
 Get the internal spectral data.
 
std::pair< BoxArray, DistributionMappinggetSpectralDataLayout () 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 *, IntVectgetSpectralData () const
 

Private Member Functions

void prepare_openbc ()
 
void backward_doit (MF &outmf, IntVect const &ngout=IntVect(0), Periodicity const &period=Periodicity::NonPeriodic(), int outcomp=0)
 
void backward_doit (cMF const &inmf, MF &outmf, IntVect const &ngout=IntVect(0), Periodicity const &period=Periodicity::NonPeriodic(), int incomp=0, int outcomp=0)
 
std::pair< Plan< T >, Plan< T > > make_c2c_plans (cMF &inout, int ndims) const
 
template<typename FA , typename RT >
std::pair< std::unique_ptr< char, DataDeleter >, std::size_t > install_raw_ptr (FA &fa, RT const *p)
 

Static Private Member Functions

static Box make_domain_x (Box const &domain)
 
static Box make_domain_y (Box const &domain)
 
static Box make_domain_z (Box const &domain)
 
static std::pair< BoxArray, DistributionMappingmake_layout_from_local_domain (std::array< int, AMREX_SPACEDIM > const &local_start, std::array< int, AMREX_SPACEDIM > const &local_size)
 

Private Attributes

Plan< T > m_fft_fwd_x {}
 
Plan< T > m_fft_bwd_x {}
 
Plan< T > m_fft_fwd_y {}
 
Plan< T > m_fft_bwd_y {}
 
Plan< T > m_fft_fwd_z {}
 
Plan< T > m_fft_bwd_z {}
 
Plan< T > m_fft_fwd_x_half {}
 
Plan< T > m_fft_bwd_x_half {}
 
std::unique_ptr< MultiBlockCommMetaDatam_cmd_x2y
 
std::unique_ptr< MultiBlockCommMetaDatam_cmd_y2x
 
std::unique_ptr< MultiBlockCommMetaDatam_cmd_y2z
 
std::unique_ptr< MultiBlockCommMetaDatam_cmd_z2y
 
std::unique_ptr< MultiBlockCommMetaDatam_cmd_x2z
 
std::unique_ptr< MultiBlockCommMetaDatam_cmd_z2x
 
std::unique_ptr< MultiBlockCommMetaDatam_cmd_x2z_half
 
std::unique_ptr< MultiBlockCommMetaDatam_cmd_z2x_half
 
Swap01 m_dtos_x2y {}
 
Swap01 m_dtos_y2x {}
 
Swap02 m_dtos_y2z {}
 
Swap02 m_dtos_z2y {}
 
RotateFwd m_dtos_x2z {}
 
RotateBwd m_dtos_z2x {}
 
MF m_rx
 
cMF m_cx
 
cMF m_cy
 
cMF m_cz
 
MF m_raw_mf
 
cMF m_raw_cmf
 
std::unique_ptr< char, DataDeleterm_data_1
 
std::unique_ptr< char, DataDeleterm_data_2
 
Box m_real_domain
 
Box m_spectral_domain_x
 
Box m_spectral_domain_y
 
Box m_spectral_domain_z
 
std::unique_ptr< R2C< T, D, C > > m_r2c_sub
 
detail::SubHelper m_sub_helper
 
Info m_info
 
bool m_do_alld_fft = false
 
bool m_slab_decomp = false
 
bool m_openbc_half = false
 

Friends

template<typename U >
class OpenBCSolver
 
template<typename U >
class Poisson
 
template<typename U >
class PoissonHybrid
 

Detailed Description

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.

Member Typedef Documentation

◆ cMF

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
using amrex::FFT::R2C< T, D, C >::cMF = FabArray<BaseFab<GpuComplex<T> > >

◆ MF

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
using amrex::FFT::R2C< T, D, C >::MF = std::conditional_t <C, cMF, std::conditional_t<std::is_same_v<T,Real>, MultiFab, FabArray<BaseFab<T> > > >

Constructor & Destructor Documentation

◆ R2C() [1/4]

template<typename T , Direction D, bool C>
amrex::FFT::R2C< T, D, C >::R2C ( Box const &  domain,
Info const &  info = Info{} 
)
explicit

Constructor.

Parameters
domainthe forward domain (i.e., the domain of the real data)
infooptional information

◆ R2C() [2/4]

template<typename T , Direction D, bool C>
amrex::FFT::R2C< T, D, C >::R2C ( std::array< int, AMREX_SPACEDIM > const &  domain_size,
Info const &  info = Info{} 
)
explicit

Constructor.

If AMREX_SPACEDIM is 3 and you want to do 2D FFT, you just need to set the size of one of the dimensions to 1.

Parameters
domain_sizesize of the forward domain (i.e., the real data domain)
infooptional information

◆ ~R2C()

template<typename T , Direction D, bool C>
amrex::FFT::R2C< T, D, C >::~R2C ( )

◆ R2C() [3/4]

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
amrex::FFT::R2C< T, D, C >::R2C ( R2C< T, D, C > const &  )
delete

◆ R2C() [4/4]

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
amrex::FFT::R2C< T, D, C >::R2C ( R2C< T, D, C > &&  )
delete

Member Function Documentation

◆ backward() [1/3]

template<typename T , Direction D, bool C>
template<Direction DIR, std::enable_if_t< DIR==Direction::backward||DIR==Direction::both, int > >
void amrex::FFT::R2C< T, D, C >::backward ( cMF const &  inmf,
MF outmf,
int  incomp = 0,
int  outcomp = 0 
)

Backward transform.

This function is not available when this class template is instantiated for forward-only transform.

Parameters
inmfinput data in FabArray<BaseFab<GpuComplex<T>>>
outmfoutput data in MultiFab or FabArray<BaseFab<float>>

◆ backward() [2/3]

template<typename T , Direction D, bool C>
template<typename CT , typename RT , Direction DIR, bool CP, std::enable_if_t<(DIR==Direction::backward||DIR==Direction::both) &&((sizeof(RT) *2==sizeof(CT) &&!CP)||(sizeof(RT)==sizeof(CT) &&CP)), int > >
void amrex::FFT::R2C< T, D, C >::backward ( CT const *  in,
RT *  out 
)

Backward transform.

This raw pointer version of backward requires setLocalDomain/getLocalDomain and setLocalSpectralDomain/getLocalSpectralDomain have been called already. Note that one is allowed to call this function multiple times after the set/get domain functions are called only once unless the domain decomposition changes. In fact, that is the preferred way because it has better performance. All processes need to call this function even if their local size is zero. If the local size is zero, one can pass nullptrs.

◆ backward() [3/3]

template<typename T , Direction D, bool C>
template<Direction DIR, std::enable_if_t< DIR==Direction::both, int > >
void amrex::FFT::R2C< T, D, C >::backward ( MF outmf,
int  outcomp = 0 
)

Backward transform.

This function is available only when this class template is instantiated for transforms in both directions.

Parameters
outmfoutput data in MultiFab or FabArray<BaseFab<float>>

◆ backward_doit() [1/2]

template<typename T , Direction D, bool C>
void amrex::FFT::R2C< T, D, C >::backward_doit ( cMF const &  inmf,
MF outmf,
IntVect const &  ngout = IntVect(0),
Periodicity const &  period = Periodicity::NonPeriodic(),
int  incomp = 0,
int  outcomp = 0 
)
private

◆ backward_doit() [2/2]

template<typename T , Direction D, bool C>
void amrex::FFT::R2C< T, D, C >::backward_doit ( MF outmf,
IntVect const &  ngout = IntVect(0),
Periodicity const &  period = Periodicity::NonPeriodic(),
int  outcomp = 0 
)
private

◆ forward() [1/3]

template<typename T , Direction D, bool C>
template<Direction DIR, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > >
void amrex::FFT::R2C< T, D, C >::forward ( MF const &  inmf,
cMF outmf,
int  incomp = 0,
int  outcomp = 0 
)

Forward transform.

This function is not available when this class template is instantiated for backward-only transform.

Parameters
inmfinput data in MultiFab or FabArray<BaseFab<float>>
outmfoutput data in FabArray<BaseFab<GpuComplex<T>>>

◆ forward() [2/3]

template<typename T , Direction D, bool C>
template<Direction DIR, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > >
void amrex::FFT::R2C< T, D, C >::forward ( MF const &  inmf,
int  incomp = 0 
)

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
inmfinput data in MultiFab or FabArray<BaseFab<float>>

◆ forward() [3/3]

template<typename T , Direction D, bool C>
template<typename RT , typename CT , Direction DIR, bool CP, std::enable_if_t<(DIR==Direction::forward||DIR==Direction::both) &&((sizeof(RT) *2==sizeof(CT) &&!CP)||(sizeof(RT)==sizeof(CT) &&CP)), int > >
void amrex::FFT::R2C< T, D, C >::forward ( RT const *  in,
CT *  out 
)

Forward transform.

This raw pointer version of forward requires setLocalDomain/getLocalDomain and setLocalSpectralDomain/getLocalSpectralDomain have been called already. Note that one is allowed to call this function multiple times after the set/get domain functions are called only once, unless the domain decomposition changes. In fact, that is the preferred way because it has better performance. All processes need to call this function even if their local size is zero. If the local size is zero, one can pass nullptrs.

◆ forwardThenBackward()

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
template<typename F , Direction DIR = D, std::enable_if_t< DIR==Direction::both, int > = 0>
void amrex::FFT::R2C< T, D, C >::forwardThenBackward ( MF const &  inmf,
MF outmf,
F const &  post_forward,
int  incomp = 0,
int  outcomp = 0 
)
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
inmfinput data in MultiFab or FabArray<BaseFab<float>>
outmfoutput data in MultiFab or FabArray<BaseFab<float>>
post_forwarda 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.

◆ getLocalDomain()

template<typename T , Direction D, bool C>
std::pair< std::array< int, AMREX_SPACEDIM >, std::array< int, AMREX_SPACEDIM > > amrex::FFT::R2C< T, D, C >::getLocalDomain ( ) const

Get local domain.

This function returns the domain decomposition chosen by AMReX. The first part of the pair is the local starting indices, and the second part is 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, setLocalDomain and getLocalDomain, should be called.

◆ getLocalSpectralDomain()

template<typename T , Direction D, bool C>
std::pair< std::array< int, AMREX_SPACEDIM >, std::array< int, AMREX_SPACEDIM > > 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.

◆ getSpectralData() [1/2]

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
template<Direction DIR = D, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > = 0>
std::pair< cMF *, IntVect > amrex::FFT::R2C< T, D, C >::getSpectralData ( ) const

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.

◆ getSpectralData() [2/2]

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
template<Direction DIR, std::enable_if_t< DIR==Direction::forward||DIR==Direction::both, int > >
std::pair< typename R2C< T, D, C >::cMF *, IntVect > amrex::FFT::R2C< T, D, C >::getSpectralData ( ) const

◆ getSpectralDataLayout()

template<typename T , Direction D, bool C>
std::pair< BoxArray, DistributionMapping > amrex::FFT::R2C< T, D, C >::getSpectralDataLayout ( ) const

Get BoxArray and DistributionMapping for spectral data.

The returned BoxArray and DistributionMapping can be used to build FabArray<BaseFab<GpuComplex<T>>> for spectral data. The returned BoxArray has the usual order of (x,y,z).

◆ install_raw_ptr()

template<typename T , Direction D, bool C>
template<typename FA , typename RT >
std::pair< std::unique_ptr< char, DataDeleter >, std::size_t > amrex::FFT::R2C< T, D, C >::install_raw_ptr ( FA &  fa,
RT const *  p 
)
private

◆ make_c2c_plans()

template<typename T , Direction D, bool C>
std::pair< Plan< T >, Plan< T > > amrex::FFT::R2C< T, D, C >::make_c2c_plans ( cMF inout,
int  ndims 
) const
private

◆ make_domain_x()

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
static Box amrex::FFT::R2C< T, D, C >::make_domain_x ( Box const &  domain)
inlinestaticprivate

◆ make_domain_y()

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
static Box amrex::FFT::R2C< T, D, C >::make_domain_y ( Box const &  domain)
inlinestaticprivate

◆ make_domain_z()

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
static Box amrex::FFT::R2C< T, D, C >::make_domain_z ( Box const &  domain)
inlinestaticprivate

◆ make_layout_from_local_domain()

template<typename T , Direction D, bool C>
std::pair< BoxArray, DistributionMapping > amrex::FFT::R2C< T, D, C >::make_layout_from_local_domain ( std::array< int, AMREX_SPACEDIM > const &  local_start,
std::array< int, AMREX_SPACEDIM > const &  local_size 
)
staticprivate

◆ operator=() [1/2]

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
R2C & amrex::FFT::R2C< T, D, C >::operator= ( R2C< T, D, C > &&  )
delete

◆ operator=() [2/2]

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
R2C & amrex::FFT::R2C< T, D, C >::operator= ( R2C< T, D, C > const &  )
delete

◆ post_forward_doit_0()

template<typename T , Direction D, bool C>
template<typename F >
void amrex::FFT::R2C< T, D, C >::post_forward_doit_0 ( F const &  post_forward)

◆ post_forward_doit_1()

template<typename T , Direction D, bool C>
template<typename F >
void amrex::FFT::R2C< T, D, C >::post_forward_doit_1 ( F const &  post_forward)

◆ prepare_openbc()

template<typename T , Direction D, bool C>
void amrex::FFT::R2C< T, D, C >::prepare_openbc ( )
private

◆ scalingFactor()

template<typename T , Direction D, bool C>
T amrex::FFT::R2C< T, D, C >::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.

◆ setLocalDomain()

template<typename T , Direction D, bool C>
void amrex::FFT::R2C< T, D, C >::setLocalDomain ( std::array< int, AMREX_SPACEDIM > const &  local_start,
std::array< int, AMREX_SPACEDIM > 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

◆ setLocalSpectralDomain()

template<typename T , Direction D, bool C>
void amrex::FFT::R2C< T, D, C >::setLocalSpectralDomain ( std::array< int, AMREX_SPACEDIM > const &  local_start,
std::array< int, AMREX_SPACEDIM > 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

Friends And Related Symbol Documentation

◆ OpenBCSolver

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
template<typename U >
friend class OpenBCSolver
friend

◆ Poisson

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
template<typename U >
friend class Poisson
friend

◆ PoissonHybrid

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
template<typename U >
friend class PoissonHybrid
friend

Member Data Documentation

◆ m_cmd_x2y

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<MultiBlockCommMetaData> amrex::FFT::R2C< T, D, C >::m_cmd_x2y
private

◆ m_cmd_x2z

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<MultiBlockCommMetaData> amrex::FFT::R2C< T, D, C >::m_cmd_x2z
private

◆ m_cmd_x2z_half

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<MultiBlockCommMetaData> amrex::FFT::R2C< T, D, C >::m_cmd_x2z_half
private

◆ m_cmd_y2x

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<MultiBlockCommMetaData> amrex::FFT::R2C< T, D, C >::m_cmd_y2x
private

◆ m_cmd_y2z

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<MultiBlockCommMetaData> amrex::FFT::R2C< T, D, C >::m_cmd_y2z
private

◆ m_cmd_z2x

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<MultiBlockCommMetaData> amrex::FFT::R2C< T, D, C >::m_cmd_z2x
private

◆ m_cmd_z2x_half

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<MultiBlockCommMetaData> amrex::FFT::R2C< T, D, C >::m_cmd_z2x_half
private

◆ m_cmd_z2y

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<MultiBlockCommMetaData> amrex::FFT::R2C< T, D, C >::m_cmd_z2y
private

◆ m_cx

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
cMF amrex::FFT::R2C< T, D, C >::m_cx
private

◆ m_cy

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
cMF amrex::FFT::R2C< T, D, C >::m_cy
private

◆ m_cz

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
cMF amrex::FFT::R2C< T, D, C >::m_cz
private

◆ m_data_1

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<char,DataDeleter> amrex::FFT::R2C< T, D, C >::m_data_1
private

◆ m_data_2

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<char,DataDeleter> amrex::FFT::R2C< T, D, C >::m_data_2
private

◆ m_do_alld_fft

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
bool amrex::FFT::R2C< T, D, C >::m_do_alld_fft = false
private

◆ m_dtos_x2y

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Swap01 amrex::FFT::R2C< T, D, C >::m_dtos_x2y {}
private

◆ m_dtos_x2z

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
RotateFwd amrex::FFT::R2C< T, D, C >::m_dtos_x2z {}
private

◆ m_dtos_y2x

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Swap01 amrex::FFT::R2C< T, D, C >::m_dtos_y2x {}
private

◆ m_dtos_y2z

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Swap02 amrex::FFT::R2C< T, D, C >::m_dtos_y2z {}
private

◆ m_dtos_z2x

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
RotateBwd amrex::FFT::R2C< T, D, C >::m_dtos_z2x {}
private

◆ m_dtos_z2y

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Swap02 amrex::FFT::R2C< T, D, C >::m_dtos_z2y {}
private

◆ m_fft_bwd_x

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Plan<T> amrex::FFT::R2C< T, D, C >::m_fft_bwd_x {}
private

◆ m_fft_bwd_x_half

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Plan<T> amrex::FFT::R2C< T, D, C >::m_fft_bwd_x_half {}
private

◆ m_fft_bwd_y

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Plan<T> amrex::FFT::R2C< T, D, C >::m_fft_bwd_y {}
private

◆ m_fft_bwd_z

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Plan<T> amrex::FFT::R2C< T, D, C >::m_fft_bwd_z {}
private

◆ m_fft_fwd_x

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Plan<T> amrex::FFT::R2C< T, D, C >::m_fft_fwd_x {}
private

◆ m_fft_fwd_x_half

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Plan<T> amrex::FFT::R2C< T, D, C >::m_fft_fwd_x_half {}
private

◆ m_fft_fwd_y

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Plan<T> amrex::FFT::R2C< T, D, C >::m_fft_fwd_y {}
private

◆ m_fft_fwd_z

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Plan<T> amrex::FFT::R2C< T, D, C >::m_fft_fwd_z {}
private

◆ m_info

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Info amrex::FFT::R2C< T, D, C >::m_info
private

◆ m_openbc_half

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
bool amrex::FFT::R2C< T, D, C >::m_openbc_half = false
private

◆ m_r2c_sub

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
std::unique_ptr<R2C<T,D,C> > amrex::FFT::R2C< T, D, C >::m_r2c_sub
private

◆ m_raw_cmf

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
cMF amrex::FFT::R2C< T, D, C >::m_raw_cmf
mutableprivate

◆ m_raw_mf

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
MF amrex::FFT::R2C< T, D, C >::m_raw_mf
mutableprivate

◆ m_real_domain

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Box amrex::FFT::R2C< T, D, C >::m_real_domain
private

◆ m_rx

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
MF amrex::FFT::R2C< T, D, C >::m_rx
private

◆ m_slab_decomp

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
bool amrex::FFT::R2C< T, D, C >::m_slab_decomp = false
private

◆ m_spectral_domain_x

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Box amrex::FFT::R2C< T, D, C >::m_spectral_domain_x
private

◆ m_spectral_domain_y

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Box amrex::FFT::R2C< T, D, C >::m_spectral_domain_y
private

◆ m_spectral_domain_z

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
Box amrex::FFT::R2C< T, D, C >::m_spectral_domain_z
private

◆ m_sub_helper

template<typename T = Real, FFT::Direction D = FFT::Direction::both, bool C = false>
detail::SubHelper amrex::FFT::R2C< T, D, C >::m_sub_helper
private

The documentation for this class was generated from the following file: