Highly Efficient FFT for Exascale: HeFFTe v2.3
heffte_c.h
1 /*
2  -- heFFTe --
3  Univ. of Tennessee, Knoxville
4  @date
5 */
6 
7 #ifndef HEFFTE_C_H
8 #define HEFFTE_C_H
9 
10 #ifdef __cplusplus
11 #error "do not use heffte_c.h with a C++ compiler"
12 #endif
13 /*
14  * C interface for HeFFTe
15  */
16 #include <mpi.h>
17 #include "heffte_c_defines.h"
18 
26 
52 int heffte_plan_create(int backend, int const inbox_low[3], int const inbox_high[3], int const *inbox_order,
53  int const outbox_low[3], int const outbox_high[3], int const *outbox_order,
54  MPI_Comm const comm, heffte_plan_options const *options, heffte_plan *plan);
55 
63 int heffte_plan_create_r2c(int backend, int const inbox_low[3], int const inbox_high[3], int const *inbox_order,
64  int const outbox_low[3], int const outbox_high[3], int const *outbox_order,
65  int r2c_direction, MPI_Comm const comm, heffte_plan_options const *options, heffte_plan *plan);
66 
76 
92 
102 int heffte_is_r2c(heffte_plan const plan);
103 
126 void heffte_forward_s2c(heffte_plan const plan, float const *input, void *output, int scale);
131 void heffte_forward_c2c(heffte_plan const plan, void const *input, void *output, int scale);
136 void heffte_forward_d2z(heffte_plan const plan, double const *input, void *output, int scale);
141 void heffte_forward_z2z(heffte_plan const plan, void const *input, void *output, int scale);
142 
164 void heffte_forward_s2c_buffered(heffte_plan const plan, float const *input, void *output, void *workspace, int scale);
169 void heffte_forward_c2c_buffered(heffte_plan const plan, void const *input, void *output, void *workspace, int scale);
174 void heffte_forward_d2z_buffered(heffte_plan const plan, double const *input, void *output, void *workspace, int scale);
179 void heffte_forward_z2z_buffered(heffte_plan const plan, void const *input, void *output, void *workspace, int scale);
180 
203 void heffte_backward_c2s(heffte_plan const plan, void const *input, float *output, int scale);
208 void heffte_backward_c2c(heffte_plan const plan, void const *input, void *output, int scale);
213 void heffte_backward_z2d(heffte_plan const plan, void const *input, double *output, int scale);
218 void heffte_backward_z2z(heffte_plan const plan, void const *input, void *output, int scale);
219 
241 void heffte_backward_c2s_buffered(heffte_plan const plan, void const *input, float *output, void *workspace, int scale);
246 void heffte_backward_c2c_buffered(heffte_plan const plan, void const *input, void *output, void *workspace, int scale);
251 void heffte_backward_z2d_buffered(heffte_plan const plan, void const *input, double *output, void *workspace, int scale);
256 void heffte_backward_z2z_buffered(heffte_plan const plan, void const *input, void *output, void *workspace, int scale);
257 
258 #endif /* HEFFTE_C_H */
void heffte_forward_z2z_buffered(heffte_plan const plan, void const *input, void *output, void *workspace, int scale)
Forward transform, buffered, double precision, complex to complex, see heffte_forward_s2c_buffered()
void heffte_backward_z2z_buffered(heffte_plan const plan, void const *input, void *output, void *workspace, int scale)
Backward transform, buffered, double precision, complex to complex, see heffte_backward_c2s_buffered(...
int heffte_is_r2c(heffte_plan const plan)
Returns 1 if heffte_plan_create_r2c() was called and 0 otherwise.
void heffte_backward_c2s_buffered(heffte_plan const plan, void const *input, float *output, void *workspace, int scale)
Wrapper around heffte::fft3d::backward() and heffte::fft3d_r2c::backward() with user allocated worksp...
void heffte_forward_c2c(heffte_plan const plan, void const *input, void *output, int scale)
Forward transform, single precision, complex to complex, see heffte_forward_s2c()
void heffte_backward_z2d_buffered(heffte_plan const plan, void const *input, double *output, void *workspace, int scale)
Backward transform, buffered, double precision, complex to real, see heffte_backward_c2s_buffered()
void heffte_forward_s2c(heffte_plan const plan, float const *input, void *output, int scale)
Wrapper around heffte::fft3d::forward() and heffte::fft3d_r2c::forward()
void heffte_backward_c2c(heffte_plan const plan, void const *input, void *output, int scale)
Backward transform, single precision, complex to complex, see heffte_backward_c2s()
void heffte_backward_z2z(heffte_plan const plan, void const *input, void *output, int scale)
Backward transform, double precision, complex to complex, see heffte_backward_c2s()
int heffte_size_outbox(heffte_plan const plan)
Wrapper around heffte::fft3d::size_outbox() and heffte::fft3d_r2c::size_outbox()
void heffte_forward_s2c_buffered(heffte_plan const plan, float const *input, void *output, void *workspace, int scale)
Wrapper around heffte::fft3d::forward() and heffte::fft3d_r2c::forward() with user allocated workspac...
int heffte_plan_create(int backend, int const inbox_low[3], int const inbox_high[3], int const *inbox_order, int const outbox_low[3], int const outbox_high[3], int const *outbox_order, MPI_Comm const comm, heffte_plan_options const *options, heffte_plan *plan)
Creates a new heffte_plan with the given backend, boxes, communicator, and options.
void heffte_backward_c2c_buffered(heffte_plan const plan, void const *input, void *output, void *workspace, int scale)
Backward transform, buffered, single precision, complex to complex, see heffte_backward_c2s_buffered(...
int heffte_get_backend(heffte_plan const plan)
Returns the backend used in the create command.
int heffte_size_inbox(heffte_plan const plan)
Wrapper around heffte::fft3d::size_inbox() and heffte::fft3d_r2c::size_inbox()
void heffte_backward_z2d(heffte_plan const plan, void const *input, double *output, int scale)
Backward transform, double precision, complex to real, see heffte_backward_c2s()
void heffte_forward_c2c_buffered(heffte_plan const plan, void const *input, void *output, void *workspace, int scale)
Forward transform, buffered, double precision, complex to complex, see heffte_forward_s2c_buffered()
int heffte_plan_create_r2c(int backend, int const inbox_low[3], int const inbox_high[3], int const *inbox_order, int const outbox_low[3], int const outbox_high[3], int const *outbox_order, int r2c_direction, MPI_Comm const comm, heffte_plan_options const *options, heffte_plan *plan)
Creates a new heffte_plan for an r2c operation with the given backend, boxes, communicator,...
void heffte_forward_d2z(heffte_plan const plan, double const *input, void *output, int scale)
Forward transform, double precision, real to complex, see heffte_forward_s2c()
int heffte_plan_destroy(heffte_plan plan)
Destory a heffte_plan, call the destructor of the internal object.
int heffte_size_workspace(heffte_plan const plan)
Wrapper around heffte::fft3d::size_workspace() and heffte::fft3d_r2c::size_workspace()
void heffte_forward_z2z(heffte_plan const plan, void const *input, void *output, int scale)
Forward transform, double precision, complex to complex, see heffte_forward_s2c()
void heffte_backward_c2s(heffte_plan const plan, void const *input, float *output, int scale)
Wrapper around heffte::fft3d::backward() and heffte::fft3d_r2c::backward()
void heffte_forward_d2z_buffered(heffte_plan const plan, double const *input, void *output, void *workspace, int scale)
Forward transform, buffered, double precision, complex to complex, see heffte_forward_s2c_buffered()
int heffte_set_default_options(int backend, heffte_plan_options *options)
Populates the fields of the options with the ones default for the backend.
scale
Indicates the scaling factor to apply on the result of an FFT operation.
Definition: heffte_fft3d.h:113
Generic wrapper around some fft3d/fft3d_r2c instance, use heffte_plan instead of this.
Definition: heffte_c_defines.h:133
Equivalent to heffte::plan_options but defined for the C API.
Definition: heffte_c_defines.h:113