Highly Efficient FFT for Exascale: HeFFTe v2.3
heffte::fft3d_r2c< backend_tag, index > Class Template Reference

Similar to heffte::fft3d, but computed fewer redundant coefficients when the input is real. More...

#include <heffte_fft3d_r2c.h>

Public Types

using backend_executor_c2c = typename one_dim_backend< backend_tag >::executor
 FFT executor for the complex-to-complex dimensions.
 
using backend_executor_r2c = typename one_dim_backend< backend_tag >::executor_r2c
 FFT executor for the real-to-complex dimension.
 
using location_tag = typename backend::buffer_traits< backend_tag >::location
 Type-tag that is either tag::cpu or tag::gpu to indicate the location of the data.
 
template<typename T >
using buffer_container = typename backend::buffer_traits< backend_tag >::template container< T >
 Alias to the container template associated with the backend (allows for RAII memory management).
 
template<typename T >
using real_buffer_container = buffer_container< typename define_standard_type< T >::type::value_type >
 Container of real values corresponding to the complex type T.
 
template<typename T >
using output_buffer_container = buffer_container< typename fft_output< T >::type >
 Container of the output type corresponding to T, see the table of compatible input and output types.
 
- Public Types inherited from heffte::backend::device_instance< backend::buffer_traits< backend_tag >::location >
using stream_type = void *
 The type for the internal stream, the cpu uses just a void pointer.
 

Public Member Functions

 fft3d_r2c (box3d< index > const inbox, box3d< index > const outbox, int r2c_direction, MPI_Comm const comm, plan_options const options=default_options< backend_tag >())
 Constructor creating a plan for FFT transform across the given communicator and using the box geometry. More...
 
 fft3d_r2c (typename backend::device_instance< location_tag >::stream_type gpu_stream, box3d< index > const inbox, box3d< index > const outbox, int r2c_direction, MPI_Comm const comm, plan_options const options=default_options< backend_tag >())
 See the documentation for fft3d::fft3d()
 
 fft3d_r2c (int il0, int il1, int il2, int ih0, int ih1, int ih2, int io0, int io1, int io2, int ol0, int ol1, int ol2, int oh0, int oh1, int oh2, int oo0, int oo1, int oo2, int r2c_direction, MPI_Comm const comm, bool use_reorder, int algorithm, bool use_pencils)
 Internal use only, used by the Fortran interface.
 
 fft3d_r2c (int il0, int il1, int il2, int ih0, int ih1, int ih2, int io0, int io1, int io2, int ol0, int ol1, int ol2, int oh0, int oh1, int oh2, int oo0, int oo1, int oo2, int r2c_direction, MPI_Comm const comm)
 Internal use only, used by the Fortran interface.
 
 fft3d_r2c (int il0, int il1, int il2, int ih0, int ih1, int ih2, int ol0, int ol1, int ol2, int oh0, int oh1, int oh2, int r2c_direction, MPI_Comm const comm)
 Internal use only, used by the Fortran interface.
 
long long size_inbox () const
 Returns the size of the inbox defined in the constructor.
 
long long size_outbox () const
 Returns the size of the outbox defined in the constructor.
 
box3d< index > inbox () const
 Returns the inbox.
 
box3d< index > outbox () const
 Returns the outbox.
 
size_t size_workspace () const
 Returns the workspace size that will be used, size is measured in complex numbers.
 
size_t size_comm_buffers () const
 Returns the size used by the communication workspace buffers (internal use).
 
template<typename input_type , typename output_type >
void forward (input_type const input[], output_type output[], scale scaling=scale::none) const
 Performs a forward Fourier transform using two arrays. More...
 
template<typename input_type , typename output_type >
void forward (input_type const input[], output_type output[], output_type workspace[], scale scaling=scale::none) const
 Overload utilizing a user provided buffer.
 
template<typename input_type , typename output_type >
void forward (int batch_size, input_type const input[], output_type output[], output_type workspace[], scale scaling=scale::none) const
 Overload utilizing a batch transform.
 
template<typename input_type , typename output_type >
void forward (int batch_size, input_type const input[], output_type output[], scale scaling=scale::none) const
 Overload utilizing a batch transform using internally allocated workspace.
 
template<typename input_type >
output_buffer_container< input_type > forward (buffer_container< input_type > const &input, scale scaling=scale::none)
 Vector variant of forward() using input and output buffer_container classes. More...
 
template<typename input_type , typename output_type >
void backward (input_type const input[], output_type output[], scale scaling=scale::none) const
 Performs a backward Fourier transform using two arrays. More...
 
template<typename input_type , typename output_type >
void backward (input_type const input[], output_type output[], input_type workspace[], scale scaling=scale::none) const
 Overload utilizing a user provided buffer.
 
template<typename input_type , typename output_type >
void backward (int batch_size, input_type const input[], output_type output[], input_type workspace[], scale scaling=scale::none) const
 Overload that performs a batch transform.
 
template<typename input_type , typename output_type >
void backward (int batch_size, input_type const input[], output_type output[], scale scaling=scale::none) const
 Overload that performs a batch transform using internally allocated workspace.
 
template<typename scalar_type >
real_buffer_container< scalar_type > backward (buffer_container< scalar_type > const &input, scale scaling=scale::none)
 Variant of backward() that uses buffer_container for RAII style of resource management.
 
double get_scale_factor (scale scaling) const
 Returns the scale factor for the given scaling.
 
- Public Member Functions inherited from heffte::backend::device_instance< backend::buffer_traits< backend_tag >::location >
 device_instance (void *=nullptr)
 Empty constructor.
 
virtual ~device_instance ()=default
 Default destructor.
 
void * stream ()
 Returns the nullptr.
 
void * stream () const
 Returns the nullptr (const case).
 
void synchronize_device () const
 Syncs the execution with the queue, no-op in the CPU case.
 

Detailed Description

template<typename backend_tag, typename index = int>
class heffte::fft3d_r2c< backend_tag, index >

Similar to heffte::fft3d, but computed fewer redundant coefficients when the input is real.

Overview
Given real input data, there is no unambiguous way to distinguish between the positive and negative direction in the complex plane; therefore, by an argument of symmetry, all complex output must come in conjugate pairs. The heffte::fft3d computes both numbers for each conjugate pair, this class aims at computing fewer redundant coefficients and thus reducing both flops and data movement. This is achieved by selecting one of the three dimensions and the data is shortened in that dimensions to contain only the unique (non-conjugate) coefficients.
Boxes and Data Distribution
Similar to heffte::fft3d the data is organized in boxes using the heffte::box3d structs; however, in the real-to-complex case the global input and output domains do not match. If the original data sits in a box {0, 0, 0}, {x, y, z}, then depending on the dimensions chosen for the shortening, the output data will form the box:
{{0, 0, 0}, {x/2 + 1, y, z}} // if chosen dimension 0
{{0, 0, 0}, {x, y/2 + 1, z}} // if chosen dimension 1
{{0, 0, 0}, {x, y, z/2 + 1}} // if chosen dimension 2
// note that x/2 indicates the standard C++ integer division
Thus, the union of the inboxes across all MPI ranks must add up to the global input box, and the union of the outboxes must add up to the shortened global box.
Compatible Types
The real-to-complex variant does not support the cases when the input is complex, the supported types are the ones with real input in the table of compatible types.

Constructor & Destructor Documentation

◆ fft3d_r2c()

template<typename backend_tag , typename index = int>
heffte::fft3d_r2c< backend_tag, index >::fft3d_r2c ( box3d< index > const  inbox,
box3d< index > const  outbox,
int  r2c_direction,
MPI_Comm const  comm,
plan_options const  options = default_options<backend_tag>() 
)
inline

Constructor creating a plan for FFT transform across the given communicator and using the box geometry.

Parameters
inboxis the box for the non-transformed data, i.e., the input for the forward() transform and the output of the backward() transform.
outboxis the box for the transformed data, i.e., the output for the forward() transform and the input of the backward() transform.
r2c_directionindicates the direction where the total set of coefficients will be reduced to hold only the non-conjugate pairs; selecting a dimension with odd number of indexes will result in (slightly) smaller final data set.
commis the MPI communicator with all ranks that will participate in the FFT.
optionsis a set of options that define the FFT plan, see heffte::plan_options for details.

Member Function Documentation

◆ forward() [1/2]

template<typename backend_tag , typename index = int>
template<typename input_type , typename output_type >
void heffte::fft3d_r2c< backend_tag, index >::forward ( input_type const  input[],
output_type  output[],
scale  scaling = scale::none 
) const
inline

Performs a forward Fourier transform using two arrays.

Template Parameters
input_typeis either float or double type.
output_typeis a type compatible with the output of a forward FFT, see the table of compatible types.
Parameters
inputis an array of size at least size_inbox() holding the input data corresponding to the inbox
outputis an array of size at least size_outbox() and will be overwritten with the result from the transform corresponding to the outbox
scalingdefines the type of scaling to apply (default no-scaling).

◆ forward() [2/2]

template<typename backend_tag , typename index = int>
template<typename input_type >
output_buffer_container<input_type> heffte::fft3d_r2c< backend_tag, index >::forward ( buffer_container< input_type > const &  input,
scale  scaling = scale::none 
)
inline

Vector variant of forward() using input and output buffer_container classes.

Returns either std::vector or heffte::cuda:vector using only the C++ standard types. Allows for more C++-ish calls and RAII memory management, see the heffte::fft3d equivalent.

Template Parameters
input_typeis either float or double.
Parameters
inputis a std::vector or heffte::cuda::vector with size at least size_inbox() corresponding to the input of forward().
scalingdefines the type of scaling to apply (default no-scaling).
Returns
std::vector or heffte::cuda::vector with entries corresponding to the output type and with size equal to size_outbox() corresponding to the output of forward().
Exceptions
std::invalid_argumentis the size of the input is less than size_inbox().

◆ backward()

template<typename backend_tag , typename index = int>
template<typename input_type , typename output_type >
void heffte::fft3d_r2c< backend_tag, index >::backward ( input_type const  input[],
output_type  output[],
scale  scaling = scale::none 
) const
inline

Performs a backward Fourier transform using two arrays.

Template Parameters
input_typeis either float or double.
output_typeis a type compatible with the output of a backward FFT, see the table of compatible types.
Parameters
inputis an array of size at least size_outbox() holding the input data corresponding to the outbox
outputis an array of size at least size_inbox() and will be overwritten with the result from the transform corresponding to the inbox
scalingdefines the type of scaling to apply (default no-scaling)

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