NFFT  3.4.1
Data Structures | Macros | Functions
NFFT - Nonequispaced fast Fourier transform

Direct and fast computation of the nonequispaced discrete Fourier transform. More...

Data Structures

struct  nfft_plan
 data structure for an NFFT (nonequispaced fast Fourier transform) plan with double precision More...
 

Macros

#define PRE_PHI_HUT   (1U<<0)
 
#define FG_PSI   (1U<<1)
 
#define PRE_LIN_PSI   (1U<<2)
 
#define PRE_FG_PSI   (1U<<3)
 
#define PRE_PSI   (1U<<4)
 
#define PRE_FULL_PSI   (1U<<5)
 
#define MALLOC_X   (1U<<6)
 
#define MALLOC_F_HAT   (1U<<7)
 
#define MALLOC_F   (1U<<8)
 
#define FFT_OUT_OF_PLACE   (1U<<9)
 
#define FFTW_INIT   (1U<<10)
 
#define PRE_ONE_PSI   (PRE_LIN_PSI| PRE_FG_PSI| PRE_PSI| PRE_FULL_PSI)
 

Functions

void nfft_trafo (nfft_plan *ths)
 
void nfft_adjoint (nfft_plan *ths)
 
void nfft_init_1d (nfft_plan *ths, int N1, int M)
 
void nfft_init_2d (nfft_plan *ths, int N1, int N2, int M)
 
void nfft_init_3d (nfft_plan *ths, int N1, int N2, int N3, int M)
 
void nfft_init (nfft_plan *ths, int d, int *N, int M)
 
void nfft_precompute_one_psi (nfft_plan *ths)
 
void nfft_precompute_full_psi (nfft_plan *ths)
 
void nfft_precompute_psi (nfft_plan *ths)
 
void nfft_precompute_lin_psi (nfft_plan *ths)
 
const char * nfft_check (nfft_plan *ths)
 
void nfft_finalize (nfft_plan *ths)
 

Detailed Description

Direct and fast computation of the nonequispaced discrete Fourier transform.

This module implements the nonequispaced fast Fourier transforms. In the following, we abbreviate the term "nonequispaced fast Fourier transform" by NFFT.

We introduce our notation and nomenclature for discrete Fourier transforms. Let the torus

\[ \mathbb{T}^d := \left\{ \mathbf{x}=\left(x_t\right)_{t=0,\dots,d-1}\in\mathbb{R}^{d}: \; - \frac{1}{2} \le x_t < \frac{1}{2},\; t=0,\dots,d-1 \right\} \]

of dimension $d$ be given. It will serve as domain from which the nonequispaced nodes $\mathbf{x}$ are taken. The sampling set is given by ${\cal X}:=\{\mathbf{x}_j \in {\mathbb T}^d: \,j=0,\dots,M-1\}$. Possible frequencies $\mathbf{k}$ are collected in the multi index set

\[ I_{\mathbf{N}} := \left\{ \mathbf{k}=\left(k_t\right)_{t=0,\dots,d-1}\in \mathbb{Z}^d: - \frac{N_t}{2} \le k_t < \frac{N_t}{2} ,\;t=0,\dots,d-1 \right\}. \]

Our concern is the computation of the nonequispaced discrete Fourier transform (NDFT)

\[ f_j = \sum_{\mathbf{k}\in I_{\mathbf{N}}} \hat{f}_{\mathbf{k}} {\rm e}^{-2\pi{\rm i} \mathbf{k}\mathbf{x}_j}, \qquad j=0,\dots,M-1. \]

The corresponding adjoint NDFT is the computation of

\[ \hat f_{\mathbf{k}}=\sum_{j=0}^{M-1} f_j {\rm e}^{+2\pi{\rm i} \mathbf{k}\mathbf{x}_j}, \qquad \mathbf{k}\in I_{\mathbf{N}}. \]

Direct implementations are given by nfft_direct_trafo and nfft_direct_adjoint taking ${\cal O}(|I_{\mathbf{N}}|M)$ floating point operations. Approximative realisations take only ${\cal O}(|I_{\mathbf{N}}|\log|I_{\mathbf{N}}|+M)$ floating point operations. These are provided by nfft_trafo and nfft_adjoint, respectively.

Macro Definition Documentation

#define PRE_PHI_HUT   (1U<<0)

If this flag is set, the deconvolution step (the multiplication with the diagonal matrix $\mathbf{D}$) uses precomputed values of the Fourier transformed window function.

See Also
nfft_init
nfft_init_advanced
nfft_init_guru
Author
Stefan Kunis

Definition at line 193 of file nfft3.h.

Referenced by construct(), fastsum_init_guru_source_nodes(), fastsum_init_guru_target_nodes(), fgt_init_guru(), glacier(), glacier_cv(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), inverse_radon_trafo(), linogram_dft(), linogram_fft(), main(), mpolar_dft(), mpolar_fft(), nfsft_init_advanced(), nfsoft_init_advanced(), nnfft_finalize(), nnfft_init(), nnfft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), and taylor_time_accuracy().

#define FG_PSI   (1U<<1)

If this flag is set, the convolution step (the multiplication with the sparse matrix $\mathbf{B}$) uses particular properties of the Gaussian window function to trade multiplications for direct calls to exponential function.

See Also
nfft_init
nfft_init_advanced
nfft_init_guru
Author
Stefan Kunis

Definition at line 194 of file nfft3.h.

Referenced by nfsoft_precompute().

#define PRE_LIN_PSI   (1U<<2)

If this flag is set, the convolution step (the multiplication with the sparse matrix $\mathbf{B}$) uses linear interpolation from a lookup table of equispaced samples of the window function instead of exact values of the window function. At the moment a table of size $(K+1)d$ is used, where $K=2^{10}(m+1)$. An estimate for the size of the lookup table with respect to the target accuracy should be implemented.

See Also
nfft_init
nfft_init_advanced
nfft_init_guru
Author
Stefan Kunis

Definition at line 195 of file nfft3.h.

Referenced by fastsum_precompute_source_nodes(), fastsum_precompute_target_nodes(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), inverse_radon_trafo(), linogram_fft(), mpolar_fft(), nnfft_finalize(), nnfft_init_guru(), polar_fft(), Radon_trafo(), and reconstruct().

#define PRE_FG_PSI   (1U<<3)

If this flag is set, the convolution step (the multiplication with the sparse matrix $\mathbf{B}$) uses particular properties of the Gaussian window function to trade multiplications for direct calls to exponential function (the remaining $2dM$ direct calls are precomputed).

See Also
nfft_init
nfft_init_advanced
nfft_init_guru
Author
Stefan Kunis

Definition at line 196 of file nfft3.h.

Referenced by taylor_time_accuracy().

#define PRE_PSI   (1U<<4)
#define PRE_FULL_PSI   (1U<<5)

If this flag is set, the convolution step (the multiplication with the sparse matrix $\mathbf{B}$) uses $(2m+2)^dM$ precomputed values of the window function, in addition indices of source and target vectors are stored.

See Also
nfft_init
nfft_init_advanced
nfft_init_guru
Author
Stefan Kunis

Definition at line 198 of file nfft3.h.

Referenced by fastsum_precompute_source_nodes(), fastsum_precompute_target_nodes(), glacier(), glacier_cv(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), inverse_radon_trafo(), linogram_fft(), mpolar_fft(), nnfft_finalize(), nnfft_init_guru(), polar_fft(), Radon_trafo(), and reconstruct().

#define MALLOC_X   (1U<<6)

If this flag is set, (de)allocation of the node vector is done.

See Also
nfft_init
nfft_init_advanced
nfft_init_guru
nfft_finalize
Author
Stefan Kunis

Definition at line 199 of file nfft3.h.

Referenced by construct(), fgt_init_guru(), glacier(), glacier_cv(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), inverse_radon_trafo(), linogram_dft(), linogram_fft(), mpolar_dft(), mpolar_fft(), nfsoft_init_advanced(), nnfft_finalize(), nnfft_init(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), and taylor_init().

#define MALLOC_F_HAT   (1U<<7)

If this flag is set, (de)allocation of the vector of Fourier coefficients is done.

See Also
nfft_init
nfft_init_advanced
nfft_init_guru
nfft_finalize
Author
Stefan Kunis

Definition at line 200 of file nfft3.h.

Referenced by construct(), fgt_init_guru(), glacier(), glacier_cv(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), inverse_radon_trafo(), linogram_dft(), linogram_fft(), mpolar_dft(), mpolar_fft(), nfsoft_init_advanced(), nnfft_finalize(), nnfft_init(), nnfft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), and taylor_init().

#define MALLOC_F   (1U<<8)

If this flag is set, (de)allocation of the vector of samples is done.

See Also
nfft_init
nfft_init_advanced
nfft_init_guru
nfft_finalize
Author
Stefan Kunis

Definition at line 201 of file nfft3.h.

Referenced by construct(), glacier(), glacier_cv(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), inverse_radon_trafo(), linogram_dft(), linogram_fft(), mpolar_dft(), mpolar_fft(), nfsoft_init_advanced(), nnfft_finalize(), nnfft_init(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), and taylor_init().

#define FFT_OUT_OF_PLACE   (1U<<9)
#define FFTW_INIT   (1U<<10)
#define PRE_ONE_PSI   (PRE_LIN_PSI| PRE_FG_PSI| PRE_PSI| PRE_FULL_PSI)

Summarises if precomputation is used within the convolution step (the multiplication with the sparse matrix $\mathbf{B}$). If testing against this flag is positive, nfft_precompute_one_psi has to be called.

See Also
nfft_init
nfft_init_advanced
nfft_init_guru
nfft_precompute_one_psi
nfft_finalize
Author
Stefan Kunis

Definition at line 206 of file nfft3.h.

Referenced by glacier(), glacier_cv(), and taylor_time_accuracy().

Function Documentation

void nfft_trafo ( nfft_plan ths)

Computes a NFFT, see the definition.

  • ths The pointer to a nfft plan
Author
Stefan Kunis, Daniel Potts

Referenced by construct(), mri_inh_2d1d_trafo(), mri_inh_3d_trafo(), nfsoft_trafo(), and nnfft_trafo().

void nfft_adjoint ( nfft_plan ths)

Computes an adjoint NFFT, see the definition.

  • ths The pointer to a nfft plan
Author
Stefan Kunis, Daniel Potts

Referenced by mri_inh_3d_adjoint(), nfsoft_adjoint(), nnfft_adjoint(), and reconstruct().

void nfft_init_1d ( nfft_plan ths,
int  N1,
int  M 
)

Initialisation of a transform plan, wrapper d=1.

  • ths The pointer to a nfft plan
  • N1 bandwidth
  • M The number of nodes
Author
Stefan Kunis, Daniel Potts
void nfft_init_2d ( nfft_plan ths,
int  N1,
int  N2,
int  M 
)

Initialisation of a transform plan, wrapper d=2.

  • ths The pointer to a nfft plan
  • N1 bandwidth
  • N2 bandwidth
  • M The number of nodes
Author
Stefan Kunis, Daniel Potts

Referenced by construct().

void nfft_init_3d ( nfft_plan ths,
int  N1,
int  N2,
int  N3,
int  M 
)

Initialisation of a transform plan, wrapper d=3.

  • ths The pointer to a nfft plan
  • N1 bandwidth
  • N2 bandwidth
  • N3 bandwidth
  • M The number of nodes
Author
Stefan Kunis, Daniel Potts
void nfft_init ( nfft_plan ths,
int  d,
int *  N,
int  M 
)

Initialisation of a transform plan, simple.

  • ths The pointer to a nfft plan
  • d The dimension
  • N The multi bandwidth
  • M The number of nodes
Author
Stefan Kunis, Daniel Potts
void nfft_precompute_one_psi ( nfft_plan ths)

Precomputation for a transform plan.

  • ths The pointer to a nfft plan
Author
Stefan Kunis

wrapper for precompute*_psi

if PRE_*_PSI is set the application program has to call this routine (after) setting the nodes x

Referenced by nfsoft_precompute().

void nfft_precompute_full_psi ( nfft_plan ths)

Superceded by nfft_precompute_one_psi.

Author
Stefan Kunis

Referenced by nnfft_precompute_full_psi(), and reconstruct().

void nfft_precompute_psi ( nfft_plan ths)

Superceded by nfft_precompute_one_psi.

Author
Stefan Kunis

Referenced by construct(), nnfft_precompute_psi(), and reconstruct().

void nfft_precompute_lin_psi ( nfft_plan ths)

Superceded by nfft_precompute_one_psi.

Author
Stefan Kunis

Referenced by nnfft_precompute_lin_psi(), and reconstruct().

void nfft_check ( nfft_plan ths)

Checks a transform plan for frequently used bad parameter.

  • ths The pointer to a nfft plan
Author
Stefan Kunis, Daniel Potts
void nfft_finalize ( nfft_plan ths)

Destroys a transform plan.

  • ths The pointer to a nfft plan
Author
Stefan Kunis, Daniel Potts

Referenced by construct(), mri_inh_2d1d_finalize(), mri_inh_3d_finalize(), nfsft_finalize(), nfsoft_finalize(), nnfft_finalize(), and reconstruct().