![]() |
|
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 ![]() | |
#define | FG_PSI (1U<< 1) |
If this flag is set, the convolution step (the multiplication with the sparse matrix ![]() | |
#define | PRE_LIN_PSI (1U<< 2) |
If this flag is set, the convolution step (the multiplication with the sparse matrix ![]() | |
#define | PRE_FG_PSI (1U<< 3) |
If this flag is set, the convolution step (the multiplication with the sparse matrix ![]() ![]() | |
#define | PRE_PSI (1U<< 4) |
If this flag is set, the convolution step (the multiplication with the sparse matrix ![]() ![]() | |
#define | PRE_FULL_PSI (1U<< 5) |
If this flag is set, the convolution step (the multiplication with the sparse matrix ![]() ![]() | |
#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 ![]() | |
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_trafo_1d (nfft_plan *ths) |
void | nfft_trafo_2d (nfft_plan *ths) |
void | nfft_trafo_3d (nfft_plan *ths) |
void | nfft_adjoint (nfft_plan *ths) |
Computes an adjoint NFFT, see the definition. | |
void | nfft_adjoint_1d (nfft_plan *ths) |
void | nfft_adjoint_2d (nfft_plan *ths) |
void | nfft_adjoint_3d (nfft_plan *ths) |
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. |
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 ndft_trafo and ndft_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 149 of file nfft3.h.
Referenced by construct(), fastsum_init_guru(), 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(), nfct_finalize(), nfsft_init_advanced(), nfsoft_init_advanced(), nfst_finalize(), 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 162 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 179 of file nfft3.h.
Referenced by fastsum_precompute(), inverse_linogram_fft(), inverse_mpolar_fft(), inverse_polar_fft(), Inverse_Radon_trafo(), linogram_fft(), mpolar_fft(), nfft_finalize(), nfft_precompute_one_psi(), nfsoft_init_guru(), 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 192 of file nfft3.h.
Referenced by nfft_finalize(), nfft_precompute_one_psi(), and 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 204 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(), nfsoft_init_advanced(), nfsoft_precompute(), nfst_finalize(), nfst_full_psi(), 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 217 of file nfft3.h.
Referenced by fastsum_precompute(), glacier(), glacier_cv(), 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(), and reconstruct().
#define MALLOC_X (1U<< 6) |
If this flag is set, (de)allocation of the node vector is done.
Definition at line 228 of file nfft3.h.
Referenced by construct(), fastsum_init_guru(), 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(), nfct_finalize(), nfft_finalize(), nfft_init(), nfsoft_init_advanced(), nfst_finalize(), 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 240 of file nfft3.h.
Referenced by construct(), fastsum_init_guru(), 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(), nfct_finalize(), nfft_finalize(), nfft_init(), nfsoft_init_advanced(), nfsoft_precompute(), nfst_finalize(), 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 251 of file nfft3.h.
Referenced by construct(), fastsum_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(), nfct_finalize(), nfft_finalize(), nfft_init(), nfsoft_init_advanced(), nfst_finalize(), 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 262 of file nfft3.h.
Referenced by construct(), fastsum_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(), nfct_finalize(), nfft_finalize(), nfft_init(), nfsft_init_advanced(), nfsoft_init_advanced(), nfst_finalize(), 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 273 of file nfft3.h.
Referenced by construct(), fastsum_init_guru(), 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(), nfct_finalize(), nfft_finalize(), nfft_init(), nfsft_init_advanced(), nfsoft_init_advanced(), nfst_finalize(), 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 288 of file nfft3.h.
Referenced by glacier(), glacier_cv(), and taylor_time_accuracy().
void ndft_trafo | ( | nfft_plan * | ths | ) |
Computes a NDFT, see the definition.
Definition at line 163 of file nfft.c.
Referenced by fgt_trafo(), linogram_dft(), mpolar_dft(), nfsft_trafo(), nfsoft_trafo(), polar_dft(), and taylor_time_accuracy().
void ndft_adjoint | ( | nfft_plan * | ths | ) |
Computes an adjoint NDFT, see the definition.
Definition at line 164 of file nfft.c.
Referenced by fgt_trafo(), nfsft_adjoint(), and nfsoft_adjoint().
void nfft_trafo | ( | nfft_plan * | ths | ) |
Computes a NFFT, see the definition.
Definition at line 2784 of file nfft.c.
References nfft_plan::d, nfft_plan::g, nfft_plan::g1, nfft_plan::g2, nfft_plan::g_hat, nfft_plan::my_fftw_plan1, and TIC.
Referenced by construct(), fastsum_trafo(), fgt_trafo(), glacier_cv(), linogram_fft(), mpolar_fft(), mri_inh_2d1d_trafo(), mri_inh_3d_trafo(), nfsoft_trafo(), nnfft_trafo(), polar_fft(), Radon_trafo(), and taylor_time_accuracy().
void nfft_adjoint | ( | nfft_plan * | ths | ) |
Computes an adjoint NFFT, see the definition.
Definition at line 2820 of file nfft.c.
References nfft_plan::d, nfft_plan::g, nfft_plan::g1, nfft_plan::g2, nfft_plan::g_hat, nfft_plan::my_fftw_plan2, and TIC.
Referenced by fastsum_trafo(), fgt_trafo(), mri_inh_2d1d_adjoint(), 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.
Definition at line 3102 of file nfft.c.
References nfft_init().
void nfft_init_2d | ( | nfft_plan * | ths, | |
int | N1, | |||
int | N2, | |||
int | M | |||
) |
Initialisation of a transform plan, wrapper d=2.
Definition at line 3110 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.
Definition at line 3119 of file nfft.c.
References nfft_init().
void nfft_init | ( | nfft_plan * | ths, | |
int | d, | |||
int * | N, | |||
int | M | |||
) |
Initialisation of a transform plan, simple.
Definition at line 3057 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_malloc(), nfft_next_power_of_2(), 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!!
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.
Definition at line 3082 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, nfft_plan::nfft_flags, and nfft_malloc().
Referenced by fastsum_init_guru(), 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(), mri_inh_2d1d_init_guru(), nfsft_init_guru(), nfsoft_init_guru(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_init(), and taylor_time_accuracy().
void nfft_precompute_one_psi | ( | nfft_plan * | ths | ) |
Precomputation for a transform plan.
if PRE_*_PSI is set the application program has to call this routine (after) setting the nodes x
Definition at line 2976 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 glacier(), glacier_cv(), nfsoft_precompute(), and taylor_time_accuracy().
void nfft_precompute_full_psi | ( | nfft_plan * | ths | ) |
Superceded by nfft_precompute_one_psi.
Definition at line 2935 of file nfft.c.
References nfft_plan::d, nfft_plan::m, 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.
Definition at line 2915 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.
Definition at line 2879 of file nfft.c.
References nfft_plan::d, 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(), nfsoft_init_guru(), 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.
Definition at line 3129 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.
Definition at line 3152 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, nfft_free(), PRE_FG_PSI, PRE_FULL_PSI, PRE_LIN_PSI, PRE_PSI, nfft_plan::psi, nfft_plan::psi_index_f, nfft_plan::psi_index_g, nfft_plan::sigma, and nfft_plan::x.
Referenced by construct(), fastsum_finalize(), fgt_finalize(), glacier(), glacier_cv(), 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(), nfsft_finalize(), nfsoft_finalize(), nnfft_finalize(), polar_dft(), polar_fft(), Radon_trafo(), reconstruct(), taylor_finalize(), and taylor_time_accuracy().