Documentation of all C++ functions¶
All functions have double-precision (finufft
) and single-precision
(finufftf
) versions. Do not forget this f
suffix in the latter case.
We group the simple and vectorized interfaces together, by each of the
nine transform types (dimensions 1,2,3, and types 1,2,3).
The guru interface functions are defined at the end.
You will also want to refer to the options
and error codes
which apply to all 46 routines.
A reminder on Fourier mode ordering;
see modeord.
For example, if N1=8
in a 1D type 1 or type 2 transform:
if
opts.modeord=0
: frequency indices are ordered-4,-3,-2,-1,0,1,2,3
(CMCL ordering)if
opts.modeord=1
: frequency indices are ordered0,1,2,3,-4,-3,-2,-1
(FFT ordering)
The orderings are related by a “fftshift”.
This holds for each dimension.
Multidimensional arrays are passed by a pointer to
a contiguous Fortran-style array, with the
“fastest” dimension x, then y (if present), then z (if present), then
transform number (if ntr>1
).
We do not use C/C++-style multidimensional arrays; this gives us the most
flexibility from several languages without loss of speed or memory
due to unnecessary array copying.
In all of the simple, vectorized, and plan functions below you may either pass NULL
as the last options
argument to use default options, or a pointer to a valid finufft_opts
struct.
In this latter case you will first need to create an options struct
then set default values by passing a pointer (here opts
) to the following:
void finufft_default_opts(finufft_opts* opts)
void finufftf_default_opts(finufft_opts* opts)
Set values in a NUFFT options struct to their default values.
Be sure to use the first version for double-precision and the second for single-precision. You may then change options with, for example, opts->debug=1;
and then pass opts
to the below routines.
Simple and vectorized interfaces¶
The “simple” interfaces (the first two listed in each block) perform
a single transform, whereas the “vectorized” (the last two listed in each block,
with the word “many” in the function name) perform ntr
transforms with the same set of nonuniform points but stacked complex strengths or coefficients vectors.
Note
The motivations for the vectorized interface (and guru interface, see below) are as follows. 1) It is more efficient to bin-sort the nonuniform points only once if there are not to change between transforms. 2) For small problems, certain start-up costs cause repeated calls to the simple interface to be slower than necessary. In particular, we note that FFTW takes around 0.1 ms per thread to look up stored wisdom, which for small problems (of order 10000 or less input and output data) can, sadly, dominate the runtime.
1D transforms¶
int finufft1d1(int64_t M, double* x, complex<double>* c, int iflag, double eps, int64_t
N1, complex<double>* f, finufft_opts* opts)
int finufftf1d1(int64_t M, float* x, complex<float>* c, int iflag, float eps, int64_t N1,
complex<float>* f, finufftf_opts* opts)
int finufft1d1many(int ntr, int64_t M, double* x, complex<double>* c, int iflag, double
eps, int64_t N1, complex<double>* f, finufft_opts* opts)
int finufftf1d1many(int ntr, int64_t M, float* x, complex<float>* c, int iflag, float
eps, int64_t N1, complex<float>* f, finufftf_opts* opts)
1D complex nonuniform FFT of type 1 (nonuniform to uniform).
Computes to precision eps, via a fast algorithm, one or more transforms of the form:
M-1
f[k1] = SUM c[j] exp(+/-i k1 x(j)) for -N1/2 <= k1 <= (N1-1)/2
j=0
Inputs:
ntr how many transforms (only for vectorized "many" functions, else ntr=1)
M number of nonuniform point sources
x nonuniform points (length M real array)
c source strengths (size M*ntr complex array)
iflag if >=0, uses +i in complex exponential, otherwise -i
eps desired relative precision; smaller is slower. This can be chosen
from 1e-1 down to ~ 1e-14 (in double precision) or 1e-6 (in single)
N1 number of output Fourier modes to be computed
opts pointer to options struct (see opts.rst), or NULL for defaults
Outputs:
f Fourier mode coefficients (size N1*ntr complex array)
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* complex arrays interleave Re, Im values, and their size is stated with
dimensions ordered fastest to slowest.
* Fourier frequency indices in each dimension i are the integers lying
in [-Ni/2, (Ni-1)/2]. See above, and modeord in opts.rst for possible orderings.
int finufft1d2(int64_t M, double* x, complex<double>* c, int iflag, double eps, int64_t
N1, complex<double>* f, finufft_opts* opts)
int finufftf1d2(int64_t M, float* x, complex<float>* c, int iflag, float eps, int64_t N1,
complex<float>* f, finufftf_opts* opts)
int finufft1d2many(int ntr, int64_t M, double* x, complex<double>* c, int iflag, double
eps, int64_t N1, complex<double>* f, finufft_opts* opts)
int finufftf1d2many(int ntr, int64_t M, float* x, complex<float>* c, int iflag, float
eps, int64_t N1, complex<float>* f, finufftf_opts* opts)
1D complex nonuniform FFT of type 2 (uniform to nonuniform).
Computes to precision eps, via a fast algorithm, one or more transforms of the form:
c[j] = SUM f[k1] exp(+/-i k1 x[j]) for j = 0,...,M-1
k1
where the sum is over integers -N1/2 <= k1 <= (N1-1)/2.
Inputs:
ntr how many transforms (only for vectorized "many" functions, else ntr=1)
M number of nonuniform point targets
x nonuniform points (length M real array)
iflag if >=0, uses +i in complex exponential, otherwise -i
eps desired relative precision; smaller is slower. This can be chosen
from 1e-1 down to ~ 1e-14 (in double precision) or 1e-6 (in single)
N1 number of input Fourier modes
f Fourier mode coefficients (size N1*ntr complex array)
opts pointer to options struct (see opts.rst), or NULL for defaults
Outputs:
c values at nonuniform point targets (size M*ntr complex array)
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* complex arrays interleave Re, Im values, and their size is stated with
dimensions ordered fastest to slowest.
* Fourier frequency indices in each dimension i are the integers lying
in [-Ni/2, (Ni-1)/2]. See above, and modeord in opts.rst for possible orderings.
int finufft1d3(int64_t M, double* x, complex<double>* c, int iflag, double eps, int64_t
N, double* s, complex<double>* f, finufft_opts* opts)
int finufftf1d3(int64_t M, float* x, complex<float>* c, int iflag, float eps, int64_t N,
float* s, complex<float>* f, finufftf_opts* opts)
int finufft1d3many(int ntr, int64_t M, double* x, complex<double>* c, int iflag, double
eps, int64_t N, double* s, complex<double>* f, finufft_opts* opts)
int finufftf1d3many(int ntr, int64_t M, float* x, complex<float>* c, int iflag, float
eps, int64_t N, float* s, complex<float>* f, finufftf_opts* opts)
1D complex nonuniform FFT of type 3 (nonuniform to nonuniform).
Computes to precision eps, via a fast algorithm, one or more transforms of the form:
M-1
f[k] = SUM c[j] exp(+-i s[k] x[j]), for k = 0,...,N-1
j=0
Inputs:
ntr how many transforms (only for vectorized "many" functions, else ntr=1)
M number of nonuniform point sources
x nonuniform points in R (length M real array)
c source strengths (size M*ntr complex array)
iflag if >=0, uses +i in complex exponential, otherwise -i
eps desired relative precision; smaller is slower. This can be chosen
from 1e-1 down to ~ 1e-14 (in double precision) or 1e-6 (in single)
N number of nonuniform frequency targets
s nonuniform frequency targets in R (length N real array)
opts pointer to options struct (see opts.rst), or NULL for defaults
Outputs:
f Fourier transform values at targets (size N*ntr complex array)
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* complex arrays interleave Re, Im values, and their size is stated with
dimensions ordered fastest to slowest.
2D transforms¶
int finufft2d1(int64_t M, double* x, double* y, complex<double>* c, int iflag, double
eps, int64_t N1, int64_t N2, complex<double>* f, finufft_opts* opts)
int finufftf2d1(int64_t M, float* x, float* y, complex<float>* c, int iflag, float eps,
int64_t N1, int64_t N2, complex<float>* f, finufftf_opts* opts)
int finufft2d1many(int ntr, int64_t M, double* x, double* y, complex<double>* c, int
iflag, double eps, int64_t N1, int64_t N2, complex<double>* f, finufft_opts* opts)
int finufftf2d1many(int ntr, int64_t M, float* x, float* y, complex<float>* c, int iflag,
float eps, int64_t N1, int64_t N2, complex<float>* f, finufftf_opts* opts)
2D complex nonuniform FFT of type 1 (nonuniform to uniform).
Computes to precision eps, via a fast algorithm, one or more transforms of the form:
M-1
f[k1,k2] = SUM c[j] exp(+/-i (k1 x[j] + k2 y[j]))
j=0
for -N1/2 <= k1 <= (N1-1)/2, -N2/2 <= k2 <= (N2-1)/2.
Inputs:
ntr how many transforms (only for vectorized "many" functions, else ntr=1)
M number of nonuniform point sources
x,y nonuniform point coordinates (length M real arrays)
c source strengths (size M*ntr complex array)
iflag if >=0, uses +i in complex exponential, otherwise -i
eps desired relative precision; smaller is slower. This can be chosen
from 1e-1 down to ~ 1e-14 (in double precision) or 1e-6 (in single)
N1 number of output Fourier modes to be computed (x direction)
N2 number of output Fourier modes to be computed (y direction)
opts pointer to options struct (see opts.rst), or NULL for defaults
Outputs:
f Fourier mode coefficients (size N1*N2*ntr complex array)
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* complex arrays interleave Re, Im values, and their size is stated with
dimensions ordered fastest to slowest.
* Fourier frequency indices in each dimension i are the integers lying
in [-Ni/2, (Ni-1)/2]. See above, and modeord in opts.rst for possible orderings.
int finufft2d2(int64_t M, double* x, double* y, complex<double>* c, int iflag, double
eps, int64_t N1, int64_t N2, complex<double>* f, finufft_opts* opts)
int finufftf2d2(int64_t M, float* x, float* y, complex<float>* c, int iflag, float eps,
int64_t N1, int64_t N2, complex<float>* f, finufftf_opts* opts)
int finufft2d2many(int ntr, int64_t M, double* x, double* y, complex<double>* c, int
iflag, double eps, int64_t N1, int64_t N2, complex<double>* f, finufft_opts* opts)
int finufftf2d2many(int ntr, int64_t M, float* x, float* y, complex<float>* c, int iflag,
float eps, int64_t N1, int64_t N2, complex<float>* f, finufftf_opts* opts)
2D complex nonuniform FFT of type 2 (uniform to nonuniform).
Computes to precision eps, via a fast algorithm, one or more transforms of the form:
c[j] = SUM f[k1,k2] exp(+/-i (k1 x[j] + k2 y[j])) for j = 0,...,M-1
k1,k2
where the sum is over integers -N1/2 <= k1 <= (N1-1)/2,
-N2/2 <= k2 <= (N2-1)/2.
Inputs:
ntr how many transforms (only for vectorized "many" functions, else ntr=1)
M number of nonuniform point targets
x,y nonuniform point coordinates (length M real arrays)
iflag if >=0, uses +i in complex exponential, otherwise -i
eps desired relative precision; smaller is slower. This can be chosen
from 1e-1 down to ~ 1e-14 (in double precision) or 1e-6 (in single)
N1 number of input Fourier modes (x direction)
N2 number of input Fourier modes (y direction)
f Fourier mode coefficients (size N1*N2*ntr complex array)
opts pointer to options struct (see opts.rst), or NULL for defaults
Outputs:
c values at nonuniform point targets (size M*ntr complex array)
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* complex arrays interleave Re, Im values, and their size is stated with
dimensions ordered fastest to slowest.
* Fourier frequency indices in each dimension i are the integers lying
in [-Ni/2, (Ni-1)/2]. See above, and modeord in opts.rst for possible orderings.
int finufft2d3(int64_t M, double* x, double* y, complex<double>* c, int iflag, double
eps, int64_t N, double* s, double* t, complex<double>* f, finufft_opts* opts)
int finufftf2d3(int64_t M, float* x, float* y, complex<float>* c, int iflag, float eps,
int64_t N, float* s, float* t, complex<float>* f, finufftf_opts* opts)
int finufft2d3many(int ntr, int64_t M, double* x, double* y, complex<double>* c, int
iflag, double eps, int64_t N, double* s, double* t, complex<double>* f, finufft_opts*
opts)
int finufftf2d3many(int ntr, int64_t M, float* x, float* y, complex<float>* c, int iflag,
float eps, int64_t N, float* s, float* t, complex<float>* f, finufftf_opts* opts)
2D complex nonuniform FFT of type 3 (nonuniform to nonuniform).
Computes to precision eps, via a fast algorithm, one or more transforms of the form:
M-1
f[k] = SUM c[j] exp(+-i (s[k] x[j] + t[k] y[j])), for k = 0,...,N-1
j=0
Inputs:
ntr how many transforms (only for vectorized "many" functions, else ntr=1)
M number of nonuniform point sources
x,y nonuniform point coordinates in R^2 (length M real arrays)
c source strengths (size M*ntr complex array)
iflag if >=0, uses +i in complex exponential, otherwise -i
eps desired relative precision; smaller is slower. This can be chosen
from 1e-1 down to ~ 1e-14 (in double precision) or 1e-6 (in single)
N number of nonuniform frequency targets
s,t nonuniform frequency target coordinates in R^2 (length N real arrays)
opts pointer to options struct (see opts.rst), or NULL for defaults
Outputs:
f Fourier transform values at targets (size N*ntr complex array)
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* complex arrays interleave Re, Im values, and their size is stated with
dimensions ordered fastest to slowest.
3D transforms¶
int finufft3d1(int64_t M, double* x, double* y, double* z, complex<double>* c, int iflag,
double eps, int64_t N1, int64_t N2, int64_t N3, complex<double>* f, finufft_opts* opts)
int finufftf3d1(int64_t M, float* x, float* y, float* z, complex<float>* c, int iflag,
float eps, int64_t N1, int64_t N2, int64_t N3, complex<float>* f, finufftf_opts* opts)
int finufft3d1many(int ntr, int64_t M, double* x, double* y, double* z, complex<double>*
c, int iflag, double eps, int64_t N1, int64_t N2, int64_t N3, complex<double>* f,
finufft_opts* opts)
int finufftf3d1many(int ntr, int64_t M, float* x, float* y, float* z, complex<float>* c,
int iflag, float eps, int64_t N1, int64_t N2, int64_t N3, complex<float>* f,
finufftf_opts* opts)
3D complex nonuniform FFT of type 1 (nonuniform to uniform).
Computes to precision eps, via a fast algorithm, one or more transforms of the form:
M-1
f[k1,k2] = SUM c[j] exp(+/-i (k1 x[j] + k2 y[j] + k3 z[j]))
j=0
for -N1/2 <= k1 <= (N1-1)/2, -N2/2 <= k2 <= (N2-1)/2, -N3/2 <= k3 <= (N3-1)/2
Inputs:
ntr how many transforms (only for vectorized "many" functions, else ntr=1)
M number of nonuniform point sources
x,y,z nonuniform point coordinates (length M real arrays)
c source strengths (size M*ntr complex array)
iflag if >=0, uses +i in complex exponential, otherwise -i
eps desired relative precision; smaller is slower. This can be chosen
from 1e-1 down to ~ 1e-14 (in double precision) or 1e-6 (in single)
N1 number of output Fourier modes to be computed (x direction)
N2 number of output Fourier modes to be computed (y direction)
N3 number of output Fourier modes to be computed (z direction)
opts pointer to options struct (see opts.rst), or NULL for defaults
Outputs:
f Fourier mode coefficients (size N1*N2*N3*ntr complex array)
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* complex arrays interleave Re, Im values, and their size is stated with
dimensions ordered fastest to slowest.
* Fourier frequency indices in each dimension i are the integers lying
in [-Ni/2, (Ni-1)/2]. See above, and modeord in opts.rst for possible orderings.
int finufft3d2(int64_t M, double* x, double* y, double* z, complex<double>* c, int iflag,
double eps, int64_t N1, int64_t N2, int64_t N3, complex<double>* f, finufft_opts* opts)
int finufftf3d2(int64_t M, float* x, float* y, float* z, complex<float>* c, int iflag,
float eps, int64_t N1, int64_t N2, int64_t N3, complex<float>* f, finufftf_opts* opts)
int finufft3d2many(int ntr, int64_t M, double* x, double* y, double* z, complex<double>*
c, int iflag, double eps, int64_t N1, int64_t N2, int64_t N3, complex<double>* f,
finufft_opts* opts)
int finufftf3d2many(int ntr, int64_t M, float* x, float* y, float* z, complex<float>* c,
int iflag, float eps, int64_t N1, int64_t N2, int64_t N3, complex<float>* f,
finufftf_opts* opts)
3D complex nonuniform FFT of type 2 (uniform to nonuniform).
Computes to precision eps, via a fast algorithm, one or more transforms of the form:
c[j] = SUM f[k1,k2,k3] exp(+/-i (k1 x[j] + k2 y[j] + k3 z[j]))
k1,k2,k3
for j = 0,...,M-1,
where the sum is over integers -N1/2 <= k1 <= (N1-1)/2,
-N2/2 <= k2 <= (N2-1)/2,
-N3/2 <= k3 <= (N3-1)/2.
Inputs:
ntr how many transforms (only for vectorized "many" functions, else ntr=1)
M number of nonuniform point targets
x,y,z nonuniform point coordinates (length M real arrays)
iflag if >=0, uses +i in complex exponential, otherwise -i
eps desired relative precision; smaller is slower. This can be chosen
from 1e-1 down to ~ 1e-14 (in double precision) or 1e-6 (in single)
N1 number of input Fourier modes (x direction)
N2 number of input Fourier modes (y direction)
N3 number of input Fourier modes (z direction)
f Fourier mode coefficients (size N1*N2*N3*ntr complex array)
opts pointer to options struct (see opts.rst), or NULL for defaults
Outputs:
c values at nonuniform point targets (size M*ntr complex array)
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* complex arrays interleave Re, Im values, and their size is stated with
dimensions ordered fastest to slowest.
* Fourier frequency indices in each dimension i are the integers lying
in [-Ni/2, (Ni-1)/2]. See above, and modeord in opts.rst for possible orderings.
int finufft3d3(int64_t M, double* x, double* y, double* z, complex<double>* c, int iflag,
double eps, int64_t N, double* s, double* t, double* u, complex<double>* f, finufft_opts*
opts)
int finufftf3d3(int64_t M, float* x, float* y, float* z, complex<float>* c, int iflag,
float eps, int64_t N, float* s, float* t, float* u, complex<float>* f, finufftf_opts*
opts)
int finufft3d3many(int ntr, int64_t M, double* x, double* y, double* z, complex<double>*
c, int iflag, double eps, int64_t N, double* s, double* t, double* u, complex<double>* f,
finufft_opts* opts)
int finufftf3d3many(int ntr, int64_t M, float* x, float* y, float* z, complex<float>* c,
int iflag, float eps, int64_t N, float* s, float* t, float* u, complex<float>* f,
finufftf_opts* opts)
3D complex nonuniform FFT of type 3 (nonuniform to nonuniform).
Computes to precision eps, via a fast algorithm, one or more transforms of the form:
M-1
f[k] = SUM c[j] exp(+-i (s[k] x[j] + t[k] y[j] + u[k] z[j])),
j=0
for k = 0,...,N-1.
Inputs:
ntr how many transforms (only for vectorized "many" functions, else ntr=1)
M number of nonuniform point sources
x,y,z nonuniform point coordinates in R^3 (length M real arrays)
c source strengths (size M*ntr complex array)
iflag if >=0, uses +i in complex exponential, otherwise -i
eps desired relative precision; smaller is slower. This can be chosen
from 1e-1 down to ~ 1e-14 (in double precision) or 1e-6 (in single)
N number of nonuniform frequency targets
s,t,u nonuniform frequency target coordinates in R^3 (length N real arrays)
opts pointer to options struct (see opts.rst), or NULL for defaults
Outputs:
f Fourier transform values at targets (size N*ntr complex array)
return value 0: success, 1: success but warning, >1: error (see error.rst)
Guru plan interface¶
This provides more flexibility than the simple or vectorized interfaces.
Any transform requires (at least)
calling the following four functions in order. However, within this
sequence one may insert repeated execute
calls, or another setpts
followed by more execute
calls, as long as the transform sizes (and number of transforms ntr
) are
consistent with those that have been set in the plan
and in setpts
.
Keep in mind that setpts
retains pointers to the user’s list of nonuniform points, rather than copying these points; thus the user must not change their nonuniform point arrays until after any execute
calls that use them.
Note
The plan
object (type finufft{f}_plan
) is an opaque pointer; the public interface specifies no more details that that. Under the hood in our library the plan happens to point to a C++ object of type finufft{f}_plan_s
, whose internal details the library user should not attempt to access, nor to rely on.
int finufft_makeplan(int type, int dim, int64_t* nmodes, int iflag, int ntr, double eps,
finufft_plan* plan, finufft_opts* opts)
int finufftf_makeplan(int type, int dim, int64_t* nmodes, int iflag, int ntr, float eps,
finufftf_plan* plan, finufftf_opts* opts)
Make a plan to perform one or more general transforms.
Under the hood, for type 1 and 2, this does FFTW planning and kernel Fourier
transform precomputation. For type 3, this does very little, since the FFT
sizes are not yet known.
Inputs:
type type of transform (1,2, or 3)
dim spatial dimension (1,2, or 3)
nmodes if type is 1 or 2, numbers of Fourier modes (length dim array),
ie, {N1} in 1D, {N1,N2} in 2D, or {N1,N2,N3} in 3D.
If type is 3, it is unused.
iflag if >=0, uses +i in complex exponential, otherwise -i
ntr how many transforms (only for vectorized "many" functions, else ntr=1)
eps desired relative precision; smaller is slower. This can be chosen
from 1e-1 down to ~ 1e-14 (in double precision) or 1e-6 (in single)
opts pointer to options struct (see opts.rst), or NULL for defaults
Outputs:
plan plan object (under the hood this is a pointer to another struct)
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* All available threads are planned by default (but see opts.nthreads)
* The vectorized (many vector) plan, ie ntrans>1, can be much faster
than repeated calls with the same nonuniform points.
Note that here the I/O data ordering is stacked rather than interleaved.
* For more details about the fields in the opts pointer, see opts.rst
int finufft_setpts(finufft_plan plan, int64_t M, double* x, double* y, double* z, int64_t
N, double* s, double* t, double* z)
int finufftf_setpts(finufftf_plan plan, int64_t M, float* x, float* y, float* z, int64_t
N, float* s, float* t, float* z)
Input nonuniform points with coordinates x (and possibly y, and possibly z),
and, if type 3, nonuniform frequency target coordinates s (and possibly t,
and possibly u), into an existing plan. If type is 1 or 2 then the last four
arguments are ignored. Unused dimensions are ignored.
Under the hood, for type 1 or 2, this routine bin-sorts the points (storing
just the permutation rather than new copies of the coordinates). For type
3 it also bin-sorts the frequencies, chooses two levels of grid sizes, then
plans the inner type 2 call (interpolation and FFTW).
Inputs:
M number of nonuniform spatial points (used by all types)
x nonuniform point x-coordinates (length M real array)
y if dim>1, nonuniform point y-coordinates (length M real array),
ignored otherwise
z if dim>2, nonuniform point z-coordinates (length M real array),
ignored otherwise
N number of nonuniform frequency targets (type 3 only, ignored
otherwise)
s nonuniform frequency x-coordinates (length N real array)
t if dim>1, nonuniform frequency y-coordinates (length N real array),
ignored otherwise
u if dim>2, nonuniform frequency z-coordinates (length N real array),
ignored otherwise
Input/Outputs:
plan plan object
Outputs:
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* The coordinates in x (and if nonempty, y and z) can be any real numbers.
For type 1 and 2 transforms, their definitions imply that that the
result is invariant to adding any multiple of 2pi to these coordinates.
(Internally, each coordinate is folded to [-pi,pi); as usual for
periodic functions, rounding errors are inevitable if much larger values
are input.) For type 1 these coordinates are "sources", whereas for type
2, they are "targets".
For type 3 the coordinates are "sources", and the "frequency targets"
s (and if nonempty, t and u) may also be any real numbers; the only
restriction for type 3 is that the product of source and target domain
sizes is not too large (it controls the internal fine grid size).
* The coordinates pointed to by any used arrays x, y, z, s, t, u must
not be changed between this call and the below execute call!
int finufft_execute(finufft_plan plan, complex<double>* c, complex<double>* f)
int finufftf_execute(finufftf_plan plan, complex<float>* c, complex<float>* f)
Perform one or more NUFFT transforms using previously entered nonuniform
points and an existing plan. To summarize, this maps
type 1: c -> f
type 2: f -> c
type 3: c -> f
Inputs:
plan plan object
Input/Outputs:
c For types 1 and 3, the input strengths at the nonuniform point
sources (size M*ntr complex array).
If type 2, the output values at the nonuniform point targets
(size M*ntr complex array).
f If type 1, the output Fourier mode coefficients (size N1*ntr or
N1*N2*ntr or N1*N2*N3*ntr complex array, when dim = 1, 2, or 3
respectively).
If type 2, the input Fourier mode coefficients (size N1*ntr or
N1*N2*ntr or N1*N2*N3*ntr complex array, when dim = 1, 2, or 3
respectively).
If type 3, the output values at the nonuniform frequency targets
(size N*ntr complex array).
Outputs:
return value 0: success, 1: success but warning, >1: error (see error.rst)
Notes:
* The contents of the arrays x, y, z, s, t, u must not have changed since
the finufft_setpts call that read them. The execution rereads them
(this way of doing business saves RAM).
* f and c are contiguous Fortran-style arrays with the transform number,
if ntr>1, being the "slowest" (outer) dimension.
int finufft_destroy(finufft_plan plan)
int finufftf_destroy(finufftf_plan plan)
Deallocate a plan object. This must be used upon clean-up, or before reusing
a plan in another call to finufft_makeplan.
Inputs/Outputs:
plan plan object
Outputs:
return value 0: success, 1: success but warning, >1: error (see error.rst)