NFFT Logo 3.0 API Reference
Main Page | Modules | Data Structures | Directories | File List | Data Fields | Globals

NFFT - Nonequispaced fast Fourier transform

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

Data Structures

struct  nfft_plan
 Structure for a NFFT plan. More...

Defines

#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.
#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.
#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.
#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).
#define PRE_PSI   (1U<< 4)
 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.
#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.
#define MALLOC_X   (1U<< 6)
 If this flag is set, (de)allocation of the node vector is done.
#define MALLOC_F_HAT   (1U<< 7)
 If this flag is set, (de)allocation of the vector of Fourier coefficients is done.
#define MALLOC_F   (1U<< 8)
 If this flag is set, (de)allocation of the vector of samples is done.
#define FFT_OUT_OF_PLACE   (1U<< 9)
 If this flag is set, FFTW uses disjoint input/output vectors.
#define FFTW_INIT   (1U<< 10)
 If this flag is set, fftw_init/fftw_finalize is called.
#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}$ ).

Functions

void ndft_trafo (nfft_plan *ths)
 Computes a NDFT, see the definition.
void ndft_adjoint (nfft_plan *ths)
 Computes an adjoint NDFT, see the definition.
void nfft_trafo (nfft_plan *ths)
 Computes a NFFT, see the definition.
void nfft_adjoint (nfft_plan *ths)
 Computes an adjoint NFFT, see the definition.
void nfft_init_1d (nfft_plan *ths, int N1, int M)
 Initialisation of a transform plan, wrapper d=1.
void nfft_init_2d (nfft_plan *ths, int N1, int N2, int M)
 Initialisation of a transform plan, wrapper d=2.
void nfft_init_3d (nfft_plan *ths, int N1, int N2, int N3, int M)
 Initialisation of a transform plan, wrapper d=3.
void nfft_init (nfft_plan *ths, int d, int *N, int M)
 Initialisation of a transform plan, simple.
void nfft_init_advanced (nfft_plan *ths, int d, int *N, int M, unsigned nfft_flags_on, unsigned nfft_flags_off)
 Initialisation of a transform plan, advanced.
void nfft_init_guru (nfft_plan *ths, int d, int *N, int M, int *n, int m, unsigned nfft_flags, unsigned fftw_flags)
 Initialisation of a transform plan, guru.
void nfft_precompute_one_psi (nfft_plan *ths)
 Precomputation for a transform plan.
void nfft_precompute_full_psi (nfft_plan *ths)
 Superceded by nfft_precompute_one_psi.
void nfft_precompute_psi (nfft_plan *ths)
 Superceded by nfft_precompute_one_psi.
void nfft_precompute_lin_psi (nfft_plan *ths)
 Superceded by nfft_precompute_one_psi.
void nfft_check (nfft_plan *ths)
 Checks a transform plan for frequently used bad parameter.
void nfft_finalize (nfft_plan *ths)
 Destroys a transform plan.

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,\hdots,d-1}\in\mathbb{R}^{d}: \; - \frac{1}{2} \le x_t < \frac{1}{2},\; t=0,\hdots,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,\hdots,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,\hdots,d-1}\in \mathbb{Z}^d: - \frac{N_t}{2} \le k_t < \frac{N_t}{2} ,\;t=0,\hdots,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,\hdots,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 ndft_trafo and ndft_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.


Define 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 92 of file nfft3.h.

Referenced by accuracy_pre_lin_psi(), construct(), fastsum_init_guru(), fgt_init_guru(), glacier(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_dft(), linogram_fft(), main(), mpolar_dft(), mpolar_fft(), nfct_finalize(), nfft_finalize(), nfft_init(), nfsft_init_advanced(), nfst_finalize(), nnfft_finalize(), nnfft_init(), nnfft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_time_accuracy(), and 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 105 of file nfft3.h.

Referenced by time_accuracy().

#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 122 of file nfft3.h.

Referenced by accuracy_pre_lin_psi(), fastsum_precompute(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_fft(), mpolar_fft(), nfft_finalize(), nfft_precompute_one_psi(), nnfft_finalize(), nnfft_init_guru(), polar_fft(), Radon_trafo(), reconstruct(), and time_accuracy().

#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 135 of file nfft3.h.

Referenced by nfft_finalize(), nfft_precompute_one_psi(), taylor_time_accuracy(), and time_accuracy().

#define PRE_PSI   (1U<< 4)
 

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.

See also:
nfft_init

nfft_init_advanced

nfft_init_guru

Author:
Stefan Kunis

Definition at line 147 of file nfft3.h.

Referenced by construct(), fastsum_init_guru(), fastsum_precompute(), fgt_init_guru(), fgt_init_node_dependent(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_dft(), linogram_fft(), main(), mpolar_dft(), mpolar_fft(), nfct_finalize(), nfft_finalize(), nfft_init(), nfft_precompute_one_psi(), nfsft_init_advanced(), nfst_finalize(), nfst_full_psi(), nnfft_finalize(), nnfft_init(), nnfft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), and time_accuracy().

#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 160 of file nfft3.h.

Referenced by fastsum_precompute(), glacier(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_fft(), mpolar_fft(), nfct_adjoint(), nfct_trafo(), nfft_finalize(), nfft_precompute_one_psi(), nfst_finalize(), nfst_precompute_psi(), nnfft_finalize(), nnfft_init_guru(), polar_fft(), Radon_trafo(), reconstruct(), and time_accuracy().

#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 171 of file nfft3.h.

Referenced by accuracy_pre_lin_psi(), construct(), fastsum_init_guru(), fgt_init_guru(), glacier(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_dft(), linogram_fft(), mpolar_dft(), mpolar_fft(), nfct_finalize(), nfft_finalize(), nfft_init(), nfst_finalize(), nnfft_finalize(), nnfft_init(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_init(), and time_accuracy().

#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 183 of file nfft3.h.

Referenced by accuracy_pre_lin_psi(), construct(), fastsum_init_guru(), fgt_init_guru(), glacier(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_dft(), linogram_fft(), mpolar_dft(), mpolar_fft(), nfct_finalize(), nfft_finalize(), nfft_init(), nfst_finalize(), nnfft_finalize(), nnfft_init(), nnfft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_init(), and time_accuracy().

#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 194 of file nfft3.h.

Referenced by accuracy_pre_lin_psi(), construct(), fastsum_init_guru(), glacier(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_dft(), linogram_fft(), mpolar_dft(), mpolar_fft(), nfct_finalize(), nfft_finalize(), nfft_init(), nfst_finalize(), nnfft_finalize(), nnfft_init(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_init(), and time_accuracy().

#define FFT_OUT_OF_PLACE   (1U<< 9)
 

If this flag is set, FFTW uses disjoint input/output vectors.

See also:
nfft_init

nfft_init_advanced

nfft_init_guru

nfft_finalize

Author:
Stefan Kunis

Definition at line 205 of file nfft3.h.

Referenced by accuracy_pre_lin_psi(), construct(), fastsum_init_guru(), glacier(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_dft(), linogram_fft(), main(), mpolar_dft(), mpolar_fft(), nfct_finalize(), nfft_finalize(), nfft_init(), nfsft_init_advanced(), nfst_finalize(), nnfft_init(), nnfft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_init(), taylor_time_accuracy(), and time_accuracy().

#define FFTW_INIT   (1U<< 10)
 

If this flag is set, fftw_init/fftw_finalize is called.

See also:
nfft_init

nfft_init_advanced

nfft_init_guru

nfft_finalize

Author:
Stefan Kunis

Definition at line 216 of file nfft3.h.

Referenced by accuracy_pre_lin_psi(), construct(), fastsum_init_guru(), fgt_init_guru(), glacier(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_dft(), linogram_fft(), main(), mpolar_dft(), mpolar_fft(), nfct_finalize(), nfft_finalize(), nfft_init(), nfsft_init_advanced(), nfst_finalize(), nnfft_init(), nnfft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_init(), taylor_time_accuracy(), and time_accuracy().

#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 231 of file nfft3.h.

Referenced by glacier(), nfsft_precompute_x(), simple_test_infft_1d(), and taylor_time_accuracy().


Function Documentation

void ndft_trafo nfft_plan ths  ) 
 

Computes a NDFT, see the definition.

< index over all nodes

< index for dimensions

< plain index for summation

< dito

< actual Fourier coefficient

< actual sample

< actual node x[d*j+t]

< multi index for summation

< sign times 2*pi*k*x

Definition at line 129 of file nfft.c.

Referenced by accuracy_pre_lin_psi(), fgt_trafo(), linogram_dft(), mpolar_dft(), ndft_time(), nfsft_trafo(), polar_dft(), taylor_time_accuracy(), and time_accuracy().

void ndft_adjoint nfft_plan ths  ) 
 

Computes an adjoint NDFT, see the definition.

< index over all nodes

< index for dimensions

< plain index for summation

< dito

< actual Fourier coefficient

< actual sample

< actual node x[d*j+t]

< multi index for summation

< sign times 2*pi*k*x

Definition at line 130 of file nfft.c.

Referenced by fgt_trafo(), and nfsft_adjoint().

void nfft_trafo nfft_plan ths  ) 
 

Computes a NFFT, see the definition.

form $ \hat g_k = \frac{\hat f_k}{c_k\left(\phi\right)} \text{ for } k \in I_N $

compute by d-variate discrete Fourier transform $ g_l = \sum_{k \in I_N} \hat g_k {\rm e}^{-2\pi {\rm i} \frac{kl}{n}} \text{ for } l \in I_n $

set $ f_j =\sum_{l \in I_n,m(x_j)} g_l \psi\left(x_j-\frac{l}{n}\right) \text{ for } j=0,\hdots,M_total-1 $

Definition at line 550 of file nfft.c.

Referenced by accuracy_pre_lin_psi(), construct(), fastsum_trafo(), fgt_trafo(), linogram_fft(), mpolar_fft(), mri_inh_2d1d_trafo(), mri_inh_3d_trafo(), ndft_time(), nfsft_trafo(), nnfft_trafo(), polar_fft(), Radon_trafo(), simple_test_infft_1d(), taylor_time_accuracy(), and time_accuracy().

void nfft_adjoint nfft_plan ths  ) 
 

Computes an adjoint NFFT, see the definition.

set $ g_l = \sum_{j=0}^{M_total-1} f_j \psi\left(x_j-\frac{l}{n}\right) \text{ for } l \in I_n,m(x_j) $

compute by d-variate discrete Fourier transform $ \hat g_k = \sum_{l \in I_n} g_l {\rm e}^{+2\pi {\rm i} \frac{kl}{n}} \text{ for } k \in I_N$

form $ \hat f_k = \frac{\hat g_k}{c_k\left(\phi\right)} \text{ for } k \in I_N $

Definition at line 579 of file nfft.c.

References nfft_plan::g, nfft_plan::g1, nfft_plan::g2, nfft_plan::g_hat, and nfft_plan::my_fftw_plan2.

Referenced by fastsum_trafo(), fgt_trafo(), mri_inh_2d1d_adjoint(), mri_inh_3d_adjoint(), nfsft_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

Definition at line 855 of file nfft.c.

References nfft_init().

Referenced by ndft_time(), and simple_test_infft_1d().

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

Definition at line 863 of file nfft.c.

References nfft_init().

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

Definition at line 872 of file nfft.c.

References nfft_init().

void nfft_init nfft_plan ths,
int  d,
int *  N,
int  M_total
 

Initialisation of a transform plan, simple.

< index over all dimensions

Definition at line 810 of file nfft.c.

References nfft_plan::d, FFT_OUT_OF_PLACE, nfft_plan::fftw_flags, FFTW_INIT, nfft_plan::M_total, MALLOC_F, MALLOC_F_HAT, MALLOC_X, nfft_plan::N, nfft_plan::n, nfft_plan::nfft_flags, nfft_next_power_of_2(), PRE_PHI_HUT, and PRE_PSI.

Referenced by nfft_init_1d(), nfft_init_2d(), and nfft_init_3d().

void nfft_init_advanced nfft_plan ths,
int  d,
int *  N,
int  M,
unsigned  nfft_flags_on,
unsigned  nfft_flags_off
 

Initialisation of a transform plan, advanced.

NOT YET IMPLEMENTED!!

  • ths The pointer to a nfft plan
  • d The dimension
  • N The multi bandwidth
  • M The number of nodes
  • nfft_flags_on NFFT flags to switch on
  • nfft_flags_off NFFT flags to switch off
Author:
Stefan Kunis, Daniel Potts

void nfft_init_guru nfft_plan ths,
int  d,
int *  N,
int  M_total,
int *  n,
int  m,
unsigned  nfft_flags,
unsigned  fftw_flags
 

Initialisation of a transform plan, guru.

< index over all dimensions

Definition at line 835 of file nfft.c.

References nfft_plan::d, nfft_plan::fftw_flags, nfft_plan::m, nfft_plan::M_total, nfft_plan::N, nfft_plan::n, and nfft_plan::nfft_flags.

Referenced by accuracy_pre_lin_psi(), construct(), fastsum_init_guru(), fgt_init_guru(), glacier(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_dft(), linogram_fft(), mpolar_dft(), mpolar_fft(), mri_inh_2d1d_init_guru(), mri_inh_3d_init_guru(), nfsft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_init(), taylor_time_accuracy(), and time_accuracy().

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

Definition at line 727 of file nfft.c.

References nfft_plan::nfft_flags, nfft_precompute_full_psi(), nfft_precompute_lin_psi(), nfft_precompute_psi(), PRE_FG_PSI, PRE_FULL_PSI, PRE_LIN_PSI, and PRE_PSI.

Referenced by accuracy_pre_lin_psi(), glacier(), nfsft_precompute_x(), simple_test_infft_1d(), taylor_time_accuracy(), and time_accuracy().

void nfft_precompute_full_psi nfft_plan ths  ) 
 

Superceded by nfft_precompute_one_psi.

< index over all dimensions

< index over all nodes

< plain index 0<=l_L<lprod

< multi index u<=l<=o

< multi index 0<=lj<u+o+1

< postfix plain index

< 'bandwidth' of matrix B

< depends on x_j

Definition at line 686 of file nfft.c.

References nfft_plan::d, nfft_plan::m, nfft_plan::M_total, nfft_plan::psi, nfft_plan::psi_index_f, and nfft_plan::psi_index_g.

Referenced by fastsum_precompute(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_fft(), mpolar_fft(), nfft_precompute_one_psi(), nnfft_precompute_full_psi(), polar_fft(), Radon_trafo(), and reconstruct().

void nfft_precompute_psi nfft_plan ths  ) 
 

Superceded by nfft_precompute_one_psi.

< index over all dimensions

< index over all nodes

< index u<=l<=o

< index 0<=lj<u+o+1

< depends on x_j

Definition at line 666 of file nfft.c.

References nfft_plan::d, nfft_plan::m, nfft_plan::M_total, nfft_plan::n, nfft_plan::psi, and nfft_plan::x.

Referenced by construct(), fastsum_precompute(), fgt_init_node_dependent(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_fft(), mpolar_fft(), nfft_precompute_one_psi(), nnfft_precompute_psi(), polar_fft(), Radon_trafo(), and reconstruct().

void nfft_precompute_lin_psi nfft_plan ths  ) 
 

Superceded by nfft_precompute_one_psi.

< index over all dimensions

< index over all nodes

< step size in [0,(m+1)/n]

Definition at line 630 of file nfft.c.

References nfft_plan::K, nfft_plan::m, nfft_plan::n, and nfft_plan::psi.

Referenced by fastsum_precompute(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_fft(), mpolar_fft(), nfft_precompute_one_psi(), nnfft_precompute_lin_psi(), polar_fft(), Radon_trafo(), 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

Definition at line 882 of file nfft.c.

References nfft_plan::d, nfft_plan::m, nfft_plan::M_total, nfft_plan::N, nfft_plan::sigma, and nfft_plan::x.

void nfft_finalize nfft_plan ths  ) 
 

Destroys a transform plan.

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

Definition at line 905 of file nfft.c.

References nfft_plan::c_phi_inv, nfft_plan::d, nfft_plan::f, nfft_plan::f_hat, FFT_OUT_OF_PLACE, FFTW_INIT, nfft_plan::g1, nfft_plan::g2, MALLOC_F, MALLOC_F_HAT, MALLOC_X, nfft_plan::my_fftw_plan1, nfft_plan::my_fftw_plan2, nfft_plan::n, nfft_plan::N, nfft_plan::nfft_flags, PRE_FG_PSI, PRE_FULL_PSI, PRE_LIN_PSI, PRE_PHI_HUT, PRE_PSI, nfft_plan::psi, nfft_plan::psi_index_f, nfft_plan::psi_index_g, nfft_plan::sigma, and nfft_plan::x.

Referenced by accuracy_pre_lin_psi(), construct(), fastsum_finalize(), fgt_finalize(), glacier(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_dft(), linogram_fft(), mpolar_dft(), mpolar_fft(), mri_inh_2d1d_finalize(), mri_inh_3d_finalize(), ndft_time(), nfsft_finalize(), nnfft_finalize(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), simple_test_infft_1d(), taylor_finalize(), taylor_time_accuracy(), and time_accuracy().


Generated on 1 Nov 2006 by Doxygen 1.4.4