FINUFFT.jl Reference

For installation and basic usage, see the README at https://github.com/ludvigak/FINUFFT.jl

For documentation of the library functions that are being called, see the FINUFFT documentation at https://finufft.readthedocs.io

Index

Types

FINUFFT.nufft_optsType
mutable struct nufft_opts    
    modeord            :: Cint
    chkbnds            :: Cint
    debug              :: Cint
    spread_debug       :: Cint
    showwarn           :: Cint
    nthreads           :: Cint
    fftw               :: Cint
    spread_sort        :: Cint
    spread_kerevalmeth :: Cint
    spread_kerpad      :: Cint
    upsampfac          :: Cdouble
    spread_thread      :: Cint
    maxbatchsize       :: Cint
    spread_nthr_atomic :: Cint
    spread_max_sp_size :: Cint
end

Options struct passed to the FINUFFT library.

Fields

This is a summary only; see FINUFFT documentation for full descriptions.

modeord :: Cint

0: CMCL-style increasing mode ordering (neg to pos), or
1: FFT-style mode ordering (affects type-1,2 only)

chkbnds :: Cint

0: don't check if input NU pts in [-3pi,3pi], 1: do

debug :: Cint

0: silent, 1: text basic timing output

spread_debug :: Cint

passed to spread_opts, 0 (no text) 1 (some) or 2 (lots)

showwarn :: Cint

Whether to print warnings to stderr. 0: silent, 1: print warnings

nthreads :: Cint

How many threads FINUFFT should use, or 0 (use max available in OMP)

fftw :: Cint

0:FFTW_ESTIMATE, or 1:FFTW_MEASURE (slow plan but faster FFTs)

spread_sort :: Cint

passed to spread_opts, 0 (don't sort) 1 (do) or 2 (heuristic)

spread_kerevalmeth :: Cint

passed to spread_opts, 0: exp(sqrt()), 1: Horner ppval (faster)

spread_kerpad :: Cint

passed to spread_opts, 0: don't pad to mult of 4, 1: do

upsampfac :: Cdouble

upsampling ratio sigma: 2.0 (standard), or 1.25 (small FFT), or
0.0 (auto).

spread_thread :: Cint

(for ntrans>1 only)
0: auto choice,
1: sequential multithreaded,
2: parallel singlethreaded spread.

maxbatchsize :: Cint

(for ntrans>1 only). max blocking size for vectorized, 0 for auto-set

spread_nthr_atomic :: Cint

if >=0, threads above which spreader OMP critical goes atomic

spread_max_sp_size :: Cint

if >0, overrides spreader (dir=1 only) max subproblem size

source

Functions

FINUFFT.finufft_destroy!Method
status = finufft_destroy!(plan::finufft_plan{T}) where T <: finufftReal

This destroys a FINUFFT plan object: it deallocates all stored FFTW plans, nonuniform point sorting arrays, kernel Fourier transforms arrays, and any other allocations, and nulls the plan pointer.

An integer status code is returned, that is 0 if successful. If one attempts to destroy an already-destroyed plan, 1 is returned (see FINUFFT documentation for finufft_destroy).

source
FINUFFT.finufft_exec!Method
finufft_exec!(plan::finufft_plan{T},
                  input::Array{Complex{T}},
                  output::Array{Complex{T}}) where T <: finufftReal

Execute single or many-vector FINUFFT transforms in a plan, with output written to preallocated array. See finufft_exec for arguments.

source
FINUFFT.finufft_execMethod
output::Array{Complex{T}} = finufft_exec(plan::finufft_plan{T},
                  input::Array{Complex{T}}) where T <: finufftReal

Execute single or many-vector FINUFFT transforms in a plan.

output = finufft_exec(plan, input)

For plan a previously created finufft_plan object also containing all needed nonuniform point coordinates, do a single (or if ntrans>1 in the plan stage, multiple) NUFFT transform(s), with the strengths or Fourier coefficient inputs vector(s) from input. The result of the transform(s) is returned as a (possibly multidimensional) array.

Inputs

- `plan`     `finufft_plan` object, already planned and containing
nonuniform points.
- `input`    strengths (types 1 or 3) or Fourier coefficients (type 2)
          vector, matrix, or array of appropriate size. For type 1 and 3,
          this is either a length-M vector (where M is the length of `xj`),
          or an `(M,ntrans)` matrix when `ntrans>1`. For type 2, in 1D this is size `(ms,)`, in 2D size `(ms,mt)`, or in 3D size `(ms,mt,mu)`, or
          each of these with an extra last dimension `ntrans` if `ntrans>1`.

Output

 `output`   vector of output strengths at targets (types 2 or 3), or array
          of Fourier coefficients (type 1), or, if `ntrans>1`, a stack of
          such vectors or arrays, of appropriate size.
          Specifically, if `ntrans=1`, for type 1, in 1D
          this is size `(ms,)`, in 2D size
          `(ms,mt)`, or in 3D size `(ms,mt,mu)`; for types 2 and 3
          it is a column vector of length `M` (the length of `xj` in type 2),
          or `nk` (the length of `s` in type 3). If `ntrans>1` it is a stack
          of such objects, ie, it has an extra last dimension `ntrans`.
source
FINUFFT.finufft_makeplanMethod
finufft_makeplan(type::Integer,
                      n_modes_or_dim::Union{Array{Int64},Integer},
                      iflag::Integer,
                      ntrans::Integer,
                      eps::Real;
                      dtype=Float64,
                      kwargs...)

Creates a finufft_plan object for the guru interface to FINUFFT, of type 1, 2 or 3, and with given numbers of Fourier modes (unless type 3).

Inputs

  • type transform type: 1, 2, or 3
  • n_modes_or_dim if type is 1 or 2, the number of Fourier modes in each dimension: ms in 1D, [ms mt] in 2D, or [ms mt mu] in 3D. Its length thus sets the dimension, which must be 1, 2 or 3. If type is 3, in contrast, its value fixes the dimension.
  • iflag if >=0, uses + sign in exponential, otherwise - sign.
  • ntrans number of transforms to compute simultaneously
  • eps relative precision requested (generally between 1e-15 and 1e-1), real, need not match type of dtype
  • dtype Float32 or Float64, plan for single precision or double precision
  • kwargs (optional): for options, see nufft_opts and https://finufft.readthedocs.io/en/latest/opts.html

Returns

  • finufft_plan struct

Examples

julia> p = finufft_makeplan(2,10,+1,1,1e-6);

creates a plan for a 1D type 2 Float64 transform with 10 Fourier modes and tolerance 1e-6.

julia> p = finufft_makeplan(1,[10 20],+1,1,1e-6);

creates a plan for a 2D type 1 Float64 transform with 10*20 Fourier modes.

julia> p = finufft_makeplan(3,2,+1,1,1e-4,dtype=Float32,nthreads=4);

creates a plan for a 2D type 3 Float32 transform with tolerance 1e-4, to use 4 threads.

source
FINUFFT.finufft_setpts!Method
finufft_setpts!(plan, xj [, yj[, zj[, s[, t[, u]]]]])

Input nonuniform points for general FINUFFT transform(s).

Given an already-planned finufft_plan, this reads in nonuniform point coordinate arrays xj (and yj if 2D or 3D, and zj if 3D), and additionally in the type 3 case, nonuniform frequency target coordinate arrays s (and t if 2D or 3D, and u if 3D). Empty arrays may be passed in the case of unused dimensions. For all types, sorting is done to internally store a reindexing of points, and for type 3 the spreading and FFTs are planned. These nonuniform points may then be used for multiple transforms.

Inputs

  • plan a finufft_plan object for one/many general nonuniform FFTs
  • xj Array{Float32} or Array{Float64}, vector of x-coords of all nonuniform points
  • yj empty (if dim<2), or vector of y-coords of all nonuniform points
  • zj empty (if dim<3), or vector of z-coords of all nonuniform points
  • s vector of x-coords of all nonuniform frequency targets
  • t empty (if dim<2), or vector of y-coords of all frequency targets
  • u empty (if dim<3), or vector of z-coords of all frequency targets
source
FINUFFT.nufft1d1!Method
nufft1d1!(xj      :: Array{Float64} or Array{Float32}, 
          cj      :: Array{ComplexF64} or Array{ComplexF32}, 
          iflag   :: Integer, 
          eps     :: Real,
          fk      :: Array{ComplexF64} or Array{ComplexF32};
          kwargs...
        )

Compute type-1 1D complex nonuniform FFT. Output written to fk. See nufft1d1.

source
FINUFFT.nufft1d1Method
nufft1d1(xj      :: Array{Float64} or Array{Float32}, 
         cj      :: Array{ComplexF64} or Array{ComplexF32}, 
         iflag   :: Integer, 
         eps     :: Real,
         ms      :: Integer;
         kwargs...
        ) -> Array{ComplexF64} or Array{ComplexF32}

Compute type-1 1D complex nonuniform FFT. This computes, to relative precision eps, via a fast algorithm:

          nj
f(k1) =  SUM c[j] exp(+/-i k1 x(j))  for -ms/2 <= k1 <= (ms-1)/2
         j=1

Inputs

  • xj locations of nonuniform sources on interval [-3pi,3pi), length nj
  • cj length-nj complex vector of source strengths. If length(cj)>nj, expects a stack of vectors (eg, a nj*ntrans matrix) each of which is transformed with the same source locations.
  • iflag if >=0, uses + sign in exponential, otherwise - sign.
  • eps relative precision requested (generally between 1e-15 and 1e-1)
  • ms number of Fourier modes computed, may be even or odd; in either case, mode range is integers lying in [-ms/2, (ms-1)/2]
  • kwargs (optional). See nufft_opts and https://finufft.readthedocs.io/en/latest/opts.html

Output

  • size (ms,) complex vector of Fourier coefficients f, or, if ntrans>1, matrix of size (ms,ntrans).
source
FINUFFT.nufft1d2!Method
nufft1d2!(xj      :: Array{Float64} or Array{Float32}, 
          cj      :: Array{ComplexF64} or Array{ComplexF32}, 
          iflag   :: Integer, 
          eps     :: Real,
          fk      :: Array{ComplexF64} or Array{ComplexF32};
          kwargs...
        )

Compute type-2 1D complex nonuniform FFT. Output written to cj. See nufft1d2.

source
FINUFFT.nufft1d2Method
nufft1d2(xj      :: Array{Float64} or Array{Float32}, 
         iflag   :: Integer, 
         eps     :: Real,
         fk      :: Array{ComplexF64} or Array{ComplexF32};
         kwargs...
        ) -> Array{ComplexF64}

Compute type-2 1D complex nonuniform FFT. This computes, to relative precision eps, via a fast algorithm:

c[j] = SUM   f[k1] exp(+/-i k1 x[j])      for j = 1,...,nj
        k1
 where sum is over -ms/2 <= k1 <= (ms-1)/2.

Input

  • xj location of nonuniform targets on interval [-3pi,3pi), length nj
  • fk complex Fourier coefficients. If a vector, length sets ms (with mode ordering given by opts.modeord). If a matrix, each column is transformed with the same nonuniform targets.
  • iflag if >=0, uses + sign in exponential, otherwise - sign.
  • eps relative precision requested (generally between 1e-15 and 1e-1)
  • kwargs (optional). See nufft_opts and https://finufft.readthedocs.io/en/latest/opts.html

Output

  • complex (nj,) vector c of answers at targets, or, if ntrans>1, matrix of size (nj,ntrans).
source
FINUFFT.nufft1d3!Method
nufft1d3!(xj      :: Array{Float64} or Array{Float32}, 
          cj      :: Array{ComplexF64} or Array{ComplexF32}, 
          iflag   :: Integer, 
          eps     :: Real,
          sk      :: Array{Float64} or Array{Float32},
          fk      :: Array{ComplexF64} or Array{ComplexF32};
          kwargs...
         )

Compute type-3 1D complex nonuniform FFT. Output written to fk. See nufft1d3.

source
FINUFFT.nufft1d3Method
nufft1d3(xj      :: Array{Float64} or Array{Float32}, 
         cj      :: Array{ComplexF64} or Array{ComplexF32}, 
         iflag   :: Integer, 
         eps     :: Real,
         sk      :: Array{Float64} or Array{Float32};
         kwargs...
        ) -> Array{ComplexF64}

Compute type-3 1D complex nonuniform FFT. This computes, to relative precision eps, via a fast algorithm:

         nj
f[k]  =  SUM   c[j] exp(+-i s[k] x[j]),      for k = 1, ..., nk
         j=1

Inputs

  • xj locations of nonuniform sources on R (real line), length-nj vector.
  • cj length-nj complex vector of source strengths. If length(cj)>nj, expects a size (nj,ntrans) matrix each column of which is transformed with the same source and target locations.
  • iflag if >=0, uses + sign in exponential, otherwise - sign.
  • eps relative precision requested (generally between 1e-15 and 1e-1)
  • sk frequency locations of nonuniform targets on R, length-nk vector.
  • kwargs (optional). See nufft_opts and https://finufft.readthedocs.io/en/latest/opts.html

Output

  • complex vector f size '(nk,)of values at targets, or, ifntrans>1, a matrix of size(nk,ntrans)`
source
FINUFFT.nufft2d1!Method
nufft2d1!(xj      :: Array{Float64} or Array{Float32}, 
          yj      :: Array{Float64} or Array{Float32}, 
          cj      :: Array{ComplexF64} or Array{ComplexF32}, 
          iflag   :: Integer, 
          eps     :: Real,
          fk      :: Array{ComplexF64} or Array{ComplexF32};
          kwargs...
        )

Compute type-1 2D complex nonuniform FFT. Output written to fk. See nufft2d1.

source
FINUFFT.nufft2d1Method
nufft2d1(xj      :: Array{Float64} or Array{Float32}
         yj      :: Array{Float64} or Array{Float32}, 
         cj      :: Array{ComplexF64} or Array{ComplexF32}, 
         iflag   :: Integer, 
         eps     :: Real,
         ms      :: Integer,
         mt      :: Integer;
         kwargs...
        ) -> Array{ComplexF64}

Compute type-1 2D complex nonuniform FFT. This computes, to relative precision eps, via a fast algorithm:

              nj
f[k1,k2] =   SUM  c[j] exp(+-i (k1 x[j] + k2 y[j]))
             j=1

for -ms/2 <= k1 <= (ms-1)/2,  -mt/2 <= k2 <= (mt-1)/2.

Inputs

  • xj,yj coordinates of nonuniform sources on the square [-3pi,3pi)^2, each a length-nj vector
  • cj length-nj complex vector of source strengths. If length(cj)>nj, expects a stack of vectors (eg, a nj*ntrans matrix) each of which is transformed with the same source locations.
  • iflag if >=0, uses + sign in exponential, otherwise - sign.
  • eps relative precision requested (generally between 1e-15 and 1e-1)
  • ms,mt number of Fourier modes requested in x & y; each may be even or odd. In either case the mode range is integers lying in [-m/2, (m-1)/2]
  • kwargs (optional), see nufft_opts and https://finufft.readthedocs.io/en/latest/opts.html

Output

  • size (ms,mt) complex matrix of Fourier coefficients f (ordering given by opts.modeord in each dimension; ms fast, mt slow), or, if ntrans>1, a array of size (ms,mt,ntrans).
source
FINUFFT.nufft2d2!Method
nufft2d2!(xj      :: Array{Float64} or Array{Float32}, 
          yj      :: Array{Float64} or Array{Float32}, 
          cj      :: Array{ComplexF64} or Array{ComplexF32}, 
          iflag   :: Integer, 
          eps     :: Real,
          fk      :: Array{ComplexF64} or Array{ComplexF32};
          kwargs...
        )

Compute type-2 2D complex nonuniform FFT. Output written to cj. See nufft2d2.

source
FINUFFT.nufft2d2Method
nufft2d2(xj      :: Array{Float64} or Array{Float32}, 
         yj      :: Array{Float64} or Array{Float32}, 
         iflag   :: Integer, 
         eps     :: Real,
         fk      :: Array{ComplexF64} or Array{ComplexF32};
         kwargs...
        ) -> Array{ComplexF64}

Compute type-2 2D complex nonuniform FFT. This computes, to relative precision eps, via a fast algorithm:

c[j] =  SUM   f[k1,k2] exp(+/-i (k1 x[j] + k2 y[j]))  for j = 1,..,nj
       k1,k2
 where sum is over -ms/2 <= k1 <= (ms-1)/2, -mt/2 <= k2 <= (mt-1)/2,

Inputs

  • xj,yj coordinates of nonuniform targets on the square [-3pi,3pi)^2, each a vector of length nj
  • fk complex Fourier coefficient matrix, whose size determines (ms,mt). (Mode ordering given by opts.modeord, in each dimension.) If a 3D array, 3rd dimension sets ntrans, and each of ntrans matrices is transformed with the same nonuniform targets.
  • iflag if >=0, uses + sign in exponential, otherwise - sign.
  • eps relative precision requested (generally between 1e-15 and 1e-1)
  • kwargs (optional). See nufft_opts and https://finufft.readthedocs.io/en/latest/opts.html

Output

  • complex size (nj,) vector c of answers at targets, or, if ntrans>1, matrix of size (nj,ntrans).
source
FINUFFT.nufft2d3!Method
nufft2d3!(xj      :: Array{Float64} or Array{Float32}, 
          yj      :: Array{Float64} or Array{Float32},
          cj      :: Array{ComplexF64} or Array{ComplexF32}, 
          iflag   :: Integer, 
          eps     :: Real,
          sk      :: Array{Float64} or Array{Float32},
          tk      :: Array{Float64} or Array{Float32},
          fk      :: Array{ComplexF64} or Array{ComplexF32};
          kwargs...
         )

Compute type-3 2D complex nonuniform FFT. Output written to 'fk'. See nufft2d3.

source
FINUFFT.nufft2d3Method
nufft2d3(xj      :: Array{Float64} or Array{Float32}, 
         yj      :: Array{Float64} or Array{Float32},
         cj      :: Array{ComplexF64} or Array{ComplexF32}, 
         iflag   :: Integer, 
         eps     :: Real,
         sk      :: Array{Float64} or Array{Float32},
         tk      :: Array{Float64} or Array{Float32};
         kwargs...
        ) -> Array{ComplexF64}

Compute type-3 2D complex nonuniform FFT. This computes, to relative precision eps, via a fast algorithm:

         nj
f[k]  =  SUM   c[j] exp(+-i (s[k] x[j] + t[k] y[j])),  for k = 1, ..., nk
         j=1

Inputs

  • xj,yj coordinates of nonuniform sources in R^2, each a length-nj vector.
  • cj complex vector (nj,) of source strengths. If length(cj)>nj, expects a (nj,ntrans) matrix, each column of which is transformed with the same source and target locations.
  • iflag if >=0, uses + sign in exponential, otherwise - sign.
  • eps relative precision requested (generally between 1e-15 and 1e-1)
  • sk,tk frequency coordinates of nonuniform targets in R^2, each a length-nk vector.
  • kwargs (optional). See nufft_opts and https://finufft.readthedocs.io/en/latest/opts.html

Output

  • complex vector size (nk,) of values at targets, or, if ntrans>1, a matrix of size (nk,ntrans)
source
FINUFFT.nufft3d1!Method
nufft3d1!(xj      :: Array{Float64} or Array{Float32}, 
          yj      :: Array{Float64} or Array{Float32}, 
          zj      :: Array{Float64} or Array{Float32}, 
          cj      :: Array{ComplexF64} or Array{ComplexF32}, 
          iflag   :: Integer, 
          eps     :: Real,
          fk      :: Array{ComplexF64} or Array{ComplexF32};
          kwargs...
        )

Compute type-1 3D complex nonuniform FFT. Output written to fk. See nufft3d1.

source
FINUFFT.nufft3d1Method
nufft3d1(xj      :: Array{Float64} or Array{Float32}, 
         yj      :: Array{Float64} or Array{Float32}, 
         zj      :: Array{Float64} or Array{Float32}, 
         cj      :: Array{ComplexF64} or Array{ComplexF32}, 
         iflag   :: Integer, 
         eps     :: Real,
         ms      :: Integer,
         mt      :: Integer,
         mu      :: Integer;
         kwargs...
        ) -> Array{ComplexF64}

Compute type-1 3D complex nonuniform FFT. This computes, to relative precision eps, via a fast algorithm:

                  nj
f[k1,k2,k3] =    SUM  c[j] exp(+-i (k1 x[j] + k2 y[j] + k3 z[j]))
                 j=1

for -ms/2 <= k1 <= (ms-1)/2,  -mt/2 <= k2 <= (mt-1)/2,
    -mu/2 <= k3 <= (mu-1)/2.

Inputs

  • xj,yj,zj coordinates of nonuniform sources on the cube [-3pi,3pi)^3, each a length-nj vector
  • cj length-nj complex vector of source strengths. If length(cj)>nj, expects a stack of vectors (eg, a nj*ntrans matrix) each of which is transformed with the same source locations.
  • iflag if >=0, uses + sign in exponential, otherwise - sign.
  • eps relative precision requested (generally between 1e-15 and 1e-1)
  • ms,mt,mu number of Fourier modes requested in x,y and z; each may be even or odd. In either case the mode range is integers lying in [-m/2, (m-1)/2]
  • kwargs (optional). See nufft_opts and https://finufft.readthedocs.io/en/latest/opts.html

Output

  • size (ms,mt,mu) complex array of Fourier coefficients f (ordering given by opts.modeord in each dimension; ms fastest, mu slowest), or, if ntrans>1, a 4D array of size (ms,mt,mu,ntrans).
source
FINUFFT.nufft3d2!Method
nufft3d2!(xj      :: Array{Float64} or Array{Float32}, 
          yj      :: Array{Float64} or Array{Float32}, 
          zj      :: Array{Float64} or Array{Float32}, 
          cj      :: Array{ComplexF64} or Array{ComplexF32}, 
          iflag   :: Integer, 
          eps     :: Real,
          fk      :: Array{ComplexF64} or Array{ComplexF32};
          kwargs...
        )

Compute type-2 3D complex nonuniform FFT. Output written to cj. See nufft3d2.

source
FINUFFT.nufft3d2Method
nufft3d2(xj      :: Array{Float64} or Array{Float32}, 
         yj      :: Array{Float64} or Array{Float32}, 
         zj      :: Array{Float64} or Array{Float32}, 
         iflag   :: Integer, 
         eps     :: Real,
         fk      :: Array{ComplexF64} or Array{ComplexF32};
         kwargs...
        ) -> Array{ComplexF64}

Compute type-2 3D complex nonuniform FFT. This computes, to relative precision eps, via a fast algorithm:

c[j] =   SUM   f[k1,k2,k3] exp(+/-i (k1 x[j] + k2 y[j] + k3 z[j]))
       k1,k2,k3
                       for j = 1,..,nj
 where sum is over -ms/2 <= k1 <= (ms-1)/2, -mt/2 <= k2 <= (mt-1)/2,
                  -mu/2 <= k3 <= (mu-1)/2.

Inputs

  • xj,yj,zj coordinates of nonuniform targets on the cube [-3pi,3pi)^3, each a vector of length nj
  • fk complex Fourier coefficient array, whose size sets (ms,mt,mu). (Mode ordering given by opts.modeord, in each dimension.) If a 4D array, 4th dimension sets ntrans, and each of ntrans 3D arrays is transformed with the same nonuniform targets.
  • iflag if >=0, uses + sign in exponential, otherwise - sign.
  • eps relative precision requested (generally between 1e-15 and 1e-1)
  • kwargs (optional). See nufft_opts and https://finufft.readthedocs.io/en/latest/opts.html

Output

  • complex vector c of size (nj,) giving answers at targets, or, if ntrans>1, matrix of size (nj,ntrans).
source
FINUFFT.nufft3d3!Method
nufft3d3!(xj      :: Array{Float64} or Array{Float32}, 
          yj      :: Array{Float64} or Array{Float32},
          zj      :: Array{Float64} or Array{Float32},
          cj      :: Array{ComplexF64} or Array{ComplexF32}, 
          iflag   :: Integer, 
          eps     :: Real,
          sk      :: Array{Float64} or Array{Float32},
          tk      :: Array{Float64} or Array{Float32},
          uk      :: Array{Float64} or Array{Float32},
          fk      :: Array{ComplexF64} or Array{ComplexF32};
          kwargs...
         )

Compute type-3 3D complex nonuniform FFT. Output written to fk. See nufft3d3.

source
FINUFFT.nufft3d3Method
nufft3d3(xj      :: Array{Float64} or Array{Float32}, 
         yj      :: Array{Float64} or Array{Float32},
         zj      :: Array{Float64} or Array{Float32},
         cj      :: Array{ComplexF64} or Array{ComplexF32}, 
         iflag   :: Integer, 
         eps     :: Real,
         sk      :: Array{Float64} or Array{Float32},
         tk      :: Array{Float64} or Array{Float32},
         uk      :: Array{Float64} or Array{Float32};
         kwargs...
        ) -> Array{ComplexF64}

Compute type-3 3D complex nonuniform FFT. This computes, to relative precision eps, via a fast algorithm:

         nj
f[k]  =  SUM   c[j] exp(+-i (s[k] x[j] + t[k] y[j] + u[k] z[j])),
         j=1
                         for k = 1, ..., nk

Inputs

  • xj,yj,zj coordinates of nonuniform sources in R^3, each a length-nj vector.
  • cj complex (nj,) vector of source strengths. If length(cj)>nj, expects a '(nj,ntrans)' matrix, each column of which is transformed with the same source and target locations.
  • iflag if >=0, uses + sign in exponential, otherwise - sign.
  • eps relative precision requested (generally between 1e-15 and 1e-1)
  • sk,tk,uk` frequency coordinates of nonuniform targets in R^3, each a length-nk vector.
  • kwargs (optional). See nufft_opts and https://finufft.readthedocs.io/en/latest/opts.html

Output

  • size (nk,) complex vector f values at targets, or, if ntrans>1, a matrix of size (nk,ntrans)
source