NFFT
3.4.1
|
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) |
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
of dimension be given. It will serve as domain from which the nonequispaced nodes are taken. The sampling set is given by . Possible frequencies are collected in the multi index set
Our concern is the computation of the nonequispaced discrete Fourier transform (NDFT)
The corresponding adjoint NDFT is the computation of
Direct implementations are given by nfft_direct_trafo and nfft_direct_adjoint taking floating point operations. Approximative realisations take only floating point operations. These are provided by nfft_trafo and nfft_adjoint, respectively.
#define PRE_PHI_HUT (1U<<0) |
If this flag is set, the deconvolution step (the multiplication with the diagonal matrix ) uses precomputed values of the Fourier transformed window function.
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 ) uses particular properties of the Gaussian window function to trade multiplications for direct calls to exponential function.
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 ) 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 is used, where . An estimate for the size of the lookup table with respect to the target accuracy should be implemented.
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 ) uses particular properties of the Gaussian window function to trade multiplications for direct calls to exponential function (the remaining direct calls are precomputed).
Definition at line 196 of file nfft3.h.
Referenced by taylor_time_accuracy().
#define PRE_PSI (1U<<4) |
If this flag is set, the convolution step (the multiplication with the sparse matrix ) uses precomputed values of the window function.
Definition at line 197 of file nfft3.h.
Referenced by construct(), fastsum_init_guru_source_nodes(), fastsum_init_guru_target_nodes(), fastsum_precompute_source_nodes(), fastsum_precompute_target_nodes(), 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(), nfsft_init_advanced(), nfsoft_init_advanced(), nfsoft_precompute(), nnfft_finalize(), nnfft_init(), nnfft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), and reconstruct().
#define PRE_FULL_PSI (1U<<5) |
If this flag is set, the convolution step (the multiplication with the sparse matrix ) uses precomputed values of the window function, in addition indices of source and target vectors are stored.
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.
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.
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.
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) |
If this flag is set, FFTW uses disjoint input/output vectors.
Definition at line 202 of file nfft3.h.
Referenced by construct(), fastsum_init_guru_source_nodes(), fastsum_init_guru_target_nodes(), 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_init(), nnfft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_init(), and taylor_time_accuracy().
#define FFTW_INIT (1U<<10) |
If this flag is set, fftw_init/fftw_finalize is called.
Definition at line 203 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_init(), nnfft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_init(), and taylor_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 ). If testing against this flag is positive, nfft_precompute_one_psi has to be called.
Definition at line 206 of file nfft3.h.
Referenced by glacier(), glacier_cv(), and taylor_time_accuracy().
void nfft_trafo | ( | nfft_plan * | ths | ) |
Computes a NFFT, see the definition.
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.
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.
void nfft_init_2d | ( | nfft_plan * | ths, |
int | N1, | ||
int | N2, | ||
int | M | ||
) |
Initialisation of a transform plan, wrapper d=2.
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.
void nfft_init | ( | nfft_plan * | ths, |
int | d, | ||
int * | N, | ||
int | M | ||
) |
Initialisation of a transform plan, simple.
void nfft_precompute_one_psi | ( | nfft_plan * | ths | ) |
Precomputation for a transform plan.
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.
Referenced by nnfft_precompute_full_psi(), and reconstruct().
void nfft_precompute_psi | ( | nfft_plan * | ths | ) |
Superceded by nfft_precompute_one_psi.
Referenced by construct(), nnfft_precompute_psi(), and reconstruct().
void nfft_precompute_lin_psi | ( | nfft_plan * | ths | ) |
Superceded by nfft_precompute_one_psi.
Referenced by nnfft_precompute_lin_psi(), and reconstruct().
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.
Referenced by construct(), mri_inh_2d1d_finalize(), mri_inh_3d_finalize(), nfsft_finalize(), nfsoft_finalize(), nnfft_finalize(), and reconstruct().