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

nfft3.h File Reference

Header file for the nfft3 library. More...

#include <complex.h>
#include <fftw3.h>

Go to the source code of this file.

Data Structures

struct  nfft_plan
 Structure for a NFFT plan. More...
struct  nfct_plan
 Structure for a transform plan. More...
struct  nfst_plan
 Structure for a transform plan. More...
struct  nnfft_plan
 Structure for a transform plan. More...
struct  nsfft_plan
 Structure for a NFFT plan. More...
struct  mri_inh_2d1d_plan
 The structure for the transform plan. More...
struct  mri_inh_3d_plan
 The structure for the transform plan. More...
struct  nfsft_plan
 Structure for a NFSFT transform plan. More...
struct  infft_plan
 Structure for an inverse transform plan. More...
struct  infct_plan
 Structure for an inverse transform plan. More...
struct  infst_plan
 Structure for an inverse transform plan. More...
struct  innfft_plan
 Structure for an inverse transform plan. More...
struct  imri_inh_2d1d_plan
 Structure for an inverse transform plan. More...
struct  imri_inh_3d_plan
 Structure for an inverse transform plan. More...
struct  infsft_plan
 Structure for an inverse transform plan. More...

Defines

#define MACRO_MV_PLAN(float_type)
 Macros for public members inherited by all plan structures. Vector of samples, \ size is M_total float types.
#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}$ ).
#define MALLOC_V   (1U<< 11)
 If this flag is set, (de)allocation of the frequency node vector is done.
#define NSDFT   (1U<< 12)
 If this flag is set, the member index_sparse_to_full is (de)allocated and initialised for the use in the routine nsdft_trafo and nsdft_adjoint.
#define NFSFT_NORMALIZED   (1U << 0)
 By default, all computations are performed with respect to the unnormalized basis functions

\[ \tilde{Y}_k^n(\vartheta,\varphi) = P_k^{|n|}(\cos\vartheta) \mathrm{e}^{\mathrm{i} n \varphi}. \]

If this flag is set, all computations are carried out using the $L_2$ - normalized basis functions

\[ Y_k^n(\vartheta,\varphi) = \sqrt{\frac{2k+1}{4\pi}} P_k^{|n|}(\cos\vartheta) \mathrm{e}^{\mathrm{i} n \varphi}. \]

.

#define NFSFT_USE_NDFT   (1U << 1)
 If this flag is set, the fast NFSFT algorithms (see nfsft_trafo, nfsft_adjoint) will use internally the exact but usually slower direct NDFT algorithm in favor of fast but approximative NFFT algorithm.
#define NFSFT_USE_DPT   (1U << 2)
 If this flag is set, the fast NFSFT algorithms (see nfsft_trafo, nfsft_adjoint) will use internally the usually slower direct DPT algorithm in favor of the fast FPT algorithm.
#define NFSFT_MALLOC_X   (1U << 3)
 If this flag is set, the init methods (see nfsft_init , nfsft_init_advanced , and nfsft_init_guru) will allocate memory and the method nfsft_finalize will free the array x for you.
#define NFSFT_MALLOC_F_HAT   (1U << 5)
 If this flag is set, the init methods (see nfsft_init , nfsft_init_advanced , and nfsft_init_guru) will allocate memory and the method nfsft_finalize will free the array f_hat for you.
#define NFSFT_MALLOC_F   (1U << 6)
 If this flag is set, the init methods (see nfsft_init , nfsft_init_advanced , and nfsft_init_guru) will allocate memory and the method nfsft_finalize will free the array f for you.
#define NFSFT_PRESERVE_F_HAT   (1U << 7)
 If this flag is set, it is guaranteed that during an execution of ndsft_trafo or nfsft_trafo the content of f_hat remains unchanged.
#define NFSFT_PRESERVE_X   (1U << 8)
 If this flag is set, it is guaranteed that during an execution of ndsft_trafo, nfsft_trafo or ndsft_adjoint, nfsft_adjoint the content of x remains unchanged.
#define NFSFT_PRESERVE_F   (1U << 9)
 If this flag is set, it is guaranteed that during an execution of ndsft_adjoint or nfsft_adjoint the content of f remains unchanged.
#define NFSFT_DESTROY_F_HAT   (1U << 10)
 If this flag is set, it is explicitely allowed that during an execution of ndsft_trafo or nfsft_trafo the content of f_hat may be changed.
#define NFSFT_DESTROY_X   (1U << 11)
 If this flag is set, it is explicitely allowed that during an execution of ndsft_trafo, nfsft_trafo or ndsft_adjoint, nfsft_adjoint the content of x may be changed.
#define NFSFT_DESTROY_F   (1U << 12)
 If this flag is set, it is explicitely allowed that during an execution of ndsft_adjoint or nfsft_adjoint the content of f may be changed.
#define NFSFT_NO_DIRECT_ALGORITHM   (1U << 13)
 If this flag is set, the transforms ndsft_trafo and ndsft_adjoint do not work.
#define NFSFT_NO_FAST_ALGORITHM   (1U << 14)
 If this flag is set, the transforms nfsft_trafo and nfsft_adjoint do not work.
#define NFSFT_ZERO_F_HAT   (1U << 16)
 If this flag is set, the transforms nfsft_adjoint and ndsft_adjoint set all unused entries in f_hat not corresponding to spherical Fourier coefficients to zero.
#define NFSFT_INDEX(k, n, plan)   ((2*(plan)->N+2)*((plan)->N-n+1)+(plan)->N+k+1)
 This helper macro expands to the index $i$ corresponding to the spherical Fourier coefficient $f_hat(k,n)$ for $0 \le k \le N$ , $-k \le n \le k$ with

\[ (N+2)(N-n+1)+N+k+1 \]

.

#define NFSFT_F_HAT_SIZE(N)   ((2*N+2)*(2*N+2))
 This helper macro expands to the logical size of a spherical Fourier coefficients array for a bandwidth N.
#define FPT_NO_FAST_ALGORITHM   (1U << 2)
 If set, TODO complete comment.
#define FPT_NO_DIRECT_ALGORITHM   (1U << 3)
 If set, TODO complete comment.
#define FPT_NO_STABILIZATION   (1U << 0)
 If set, no stabilization will be used.
#define FPT_PERSISTENT_DATA   (1U << 4)
 If set, TODO complete comment.
#define FPT_FUNCTION_VALUES   (1U << 5)
 If set, the output are function values at Chebyshev nodes rather than Chebyshev coefficients.
#define FPT_AL_SYMMETRY   (1U << 6)
 TODO Don't use this flag!
#define LANDWEBER   (1U<< 0)
 If this flag is set, the Landweber (Richardson) iteration is used to compute an inverse transform.
#define STEEPEST_DESCENT   (1U<< 1)
 If this flag is set, the method of steepest descent (gradient) is used to compute an inverse transform.
#define CGNR   (1U<< 2)
 If this flag is set, the conjugate gradient method for the normal equation of first kind is used to compute an inverse transform.
#define CGNE   (1U<< 3)
 If this flag is set, the conjugate gradient method for the normal equation of second kind is used to compute an inverse transform.
#define NORMS_FOR_LANDWEBER   (1U<< 4)
 If this flag is set, the Landweber iteration updates the member dot_r_iter.
#define PRECOMPUTE_WEIGHT   (1U<< 5)
 If this flag is set, the samples are weighted, eg to cope with varying sampling density.
#define PRECOMPUTE_DAMP   (1U<< 6)
 If this flag is set, the Fourier coefficients are damped, eg to favour fast decaying coefficients.
#define MACRO_SOLVER_PLAN(MV, FLT, FLT_TYPE)
 Complete macro for mangling an inverse transform.

Typedefs

typedef fpt_set_s_fpt_set
 A set of precomputed data for a set of DPT transforms of equal maximum length.

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.
void nfct_init_1d (nfct_plan *ths_plan, int N0, int M_total)
 Creates a 1-dimensional transform plan.
void nfct_init_2d (nfct_plan *ths_plan, int N0, int N1, int M_total)
 Creates a 3-dimensional transform plan.
void nfct_init_3d (nfct_plan *ths_plan, int N0, int N1, int N2, int M_total)
 Creates a 3-dimensional transform plan.
void nfct_init (nfct_plan *ths_plan, int d, int *N, int M_total)
 Creates a d-dimensional transform plan.
void nfct_init_guru (nfct_plan *ths_plan, int d, int *N, int M_total, int *n, int m, unsigned nfct_flags, unsigned fftw_flags)
 Creates a d-dimensional transform plan.
void nfct_precompute_psi (nfct_plan *ths_plan)
 precomputes the values psi if the PRE_PSI is set the application program has to call this routine after setting the nodes this_plan->x
void nfct_trafo (nfct_plan *ths_plan)
 executes a NFCT (approximate,fast), computes for $j=0,...,M\_total-1$ $f_j^C(x_j) = \sum_{k \in I_0^{N,d}} \hat{f}_k^C * cos(2 \pi k x_j)$
void ndct_trafo (nfct_plan *ths_plan)
 executes a NDCT (exact,slow), computes for $j=0,...,M\_total-1$ $f_j^C(x_j) = \sum_{k \in I_0^{N,d}} \hat{f}_k^C * cos(2 \pi k x_j)$
void nfct_adjoint (nfct_plan *ths_plan)
 executes a transposed NFCT (approximate,fast), computes for $k \in I_0^{N,d}$ $h^C(k) = \sum_{j \in I_0^{(M\_total,1)}} f_j^C * cos(2 \pi k x_j)$
void ndct_adjoint (nfct_plan *ths_plan)
 executes a direct transposed NDCT (exact,slow), computes for $k \in I_0^{N,d}$ $h^C(k) = \sum_{j \in I_0^{(M\_total,1)}} f_j^C * cos(2 \pi k x_j)$
void nfct_finalize (nfct_plan *ths_plan)
 Destroys a plan.
double nfct_phi_hut (nfct_plan *ths_plan, int k, int d)
 do some adjustments (N,n) then compute PHI_HUT
double nfct_phi (nfct_plan *ths_plan, double x, int d)
 do some adjustments (N,n) then compute PHI
int nfct_fftw_2N (int n)
 returns 2(n-1), fftw related issue
int nfct_fftw_2N_rev (int n)
 returns 0.5n+1, fftw related issue
void nfst_init_1d (nfst_plan *ths_plan, int N0, int M_total)
 Creates a 1-dimensional transform plan.
void nfst_init_2d (nfst_plan *ths_plan, int N0, int N1, int M_total)
 Creates a 3-dimensional transform plan.
void nfst_init_3d (nfst_plan *ths_plan, int N0, int N1, int N2, int M_total)
 Creates a 3-dimensional transform plan.
void nfst_init (nfst_plan *ths_plan, int d, int *N, int M_total)
 Creates a d-dimensional transform plan.
void nfst_init_m (nfst_plan *ths_plan, int d, int *N, int M_total, int m)
 Creates a d-dimensional transform plan with pcific m.
void nfst_init_guru (nfst_plan *ths_plan, int d, int *N, int M_total, int *n, int m, unsigned nfst_flags, unsigned fftw_flags)
 Creates a d-dimensional transform plan.
void nfst_precompute_psi (nfst_plan *ths_plan)
 precomputes the values psi if the PRE_PSI is set the application program has to call this routine after setting the nodes this_plan->x
void nfst_trafo (nfst_plan *ths_plan)
 executes a NFST (approximate,fast), computes for $j=0,...,M\_total-1$ $f_j^S(x_j) = \sum_{k \in I_1^{N,d}} \hat{f}_k^S * sin(2 \pi k x_j)$
void ndst_trafo (nfst_plan *ths_plan)
 executes a NDST (exact,slow), computes for $j=0,...,M\_total-1$ $f_j^S(x_j) = \sum_{k \in I_1^{N,d}} \hat{f}_k^S * sin(2 \pi k x_j)$
void nfst_adjoint (nfst_plan *ths_plan)
 executes a transposed NFST (approximate,fast), computes for $k \in I_1^{N,d}$ $h^S(k) = \sum_{j \in I_0^{M\_total,1}} f_j^S * cos(2 \pi k x_j)$
void ndst_adjoint (nfst_plan *ths_plan)
 executes a direct transposed NDST (exact,slow), computes for $k \in I_1^{N,d}$ $h^S(k) = \sum_{j \in I_0^{M\_total,1}} f_j^S * cos(2 \pi k x_j)$
void nfst_finalize (nfst_plan *ths_plan)
 Destroys a plan.
void nfst_full_psi (nfst_plan *ths_plan, double eps)
 more memory usage, a bit faster
double nfst_phi_hut (nfst_plan *ths_plan, int k, int d)
 do some adjustments (N,n) then compute PHI_HUT
double nfst_phi (nfst_plan *ths_plan, double x, int d)
 do some adjustments (N,n) then compute PHI
int nfst_fftw_2N (int n)
 returns 2(n+1), fftw related issue
int nfst_fftw_2N_rev (int n)
 returns 0.5n-1, fftw related issue
void nnfft_init (nnfft_plan *ths_plan, int d, int N_total, int M_total, int *N)
 Creates a transform plan.
void nnfft_init_guru (nnfft_plan *ths_plan, int d, int N_total, int M_total, int *N, int *N1, int m, unsigned nnfft_flags)
 Creates a transform plan.
void nndft_trafo (nnfft_plan *ths_plan)
 Executes a direct NNDFT, i.e.
void nndft_adjoint (nnfft_plan *ths_plan)
 Executes a direct adjoint NNDFT, i.e.
void nnfft_trafo (nnfft_plan *ths_plan)
 Executes a NNFFT, i.e.
void nnfft_adjoint (nnfft_plan *ths_plan)
 Executes a adjoint NNFFT, i.e.
void nnfft_precompute_lin_psi (nnfft_plan *ths_plan)
 Precomputation for a transform plan.
void nnfft_precompute_psi (nnfft_plan *ths_plan)
 Precomputation for a transform plan.
void nnfft_precompute_full_psi (nnfft_plan *ths_plan)
 Precomputation for a transform plan.
void nnfft_precompute_phi_hut (nnfft_plan *ths_plan)
 Precomputation for a transform plan.
void nnfft_finalize (nnfft_plan *ths_plan)
 Destroys a plan.
void nsdft_trafo (nsfft_plan *ths)
 Executes an NSDFT, computes for $j=0,\hdots,M-1$ :

\[ f_j = \sum_{k\in H_N^d}\hat f_k {\rm e}^{-2\pi{\rm\scriptsize i}k x_j} \]

.

void nsdft_adjoint (nsfft_plan *ths)
 Executes an adjoint NSFFT, computes for $k\in H_N^d$ :

\[ \hat f_k = \sum_{j=0,\hdots,M-1} f_j {\rm e}^{+2\pi{\rm\scriptsize i}k x_j} \]

.

void nsfft_trafo (nsfft_plan *ths)
 Executes an NSFFT, computes fast and approximate for $j=0,\hdots,M-1$ :

\[ f_j = \sum_{k\in H_N^d}\hat f_k {\rm e}^{-2\pi{\rm\scriptsize i}k x_j} \]

.

void nsfft_adjoint (nsfft_plan *ths)
 Executes an adjoint NSFFT, computes fast and approximate for $k\in H_N^d$ :

\[ \hat f_k = \sum_{j=0,\hdots,M-1} f_j {\rm e}^{+2\pi{\rm\scriptsize i}k x_j} \]

.

void nsfft_cp (nsfft_plan *ths, nfft_plan *ths_nfft)
 Copy coefficients from nsfft plan to a nfft plan.
void nsfft_init_random_nodes_coeffs (nsfft_plan *ths)
 Initialisation of pseudo random nodes and coefficients.
void nsfft_init (nsfft_plan *ths, int d, int J, int M, int m, unsigned flags)
 Initialisation of a transform plan.
void nsfft_finalize (nsfft_plan *ths)
 Destroys a transform plan.
void mri_inh_2d1d_trafo (mri_inh_2d1d_plan *ths)
 Executes a mri transformation considering the field inhomogeneity with the 2d1d method, i.e.
void mri_inh_2d1d_adjoint (mri_inh_2d1d_plan *ths)
 Executes an adjoint mri transformation considering the field inhomogeneity with the 2d1d method, i.e.
void mri_inh_2d1d_init_guru (mri_inh_2d1d_plan *ths, int *N, int M, int *n, int m, double sigma, unsigned nfft_flags, unsigned fftw_flags)
 Creates a transform plan.
void mri_inh_2d1d_finalize (mri_inh_2d1d_plan *ths)
 Destroys a plan.
void mri_inh_3d_trafo (mri_inh_3d_plan *ths)
 Executes a mri transformation considering the field inhomogeneity with the 3d method, i.e.
void mri_inh_3d_adjoint (mri_inh_3d_plan *ths)
 Executes an adjoint mri transformation considering the field inhomogeneity with the 3d method, i.e.
void mri_inh_3d_init_guru (mri_inh_3d_plan *ths, int *N, int M, int *n, int m, double sigma, unsigned nfft_flags, unsigned fftw_flags)
void mri_inh_3d_finalize (mri_inh_3d_plan *ths)
 Destroys a plan.
void nfsft_init (nfsft_plan *plan, int N, int M)
 Creates a transform plan.
void nfsft_init_advanced (nfsft_plan *plan, int N, int M, unsigned int nfsft_flags)
 Creates a transform plan.
void nfsft_init_guru (nfsft_plan *plan, int N, int M, unsigned int nfsft_flags, int nfft_flags, int nfft_cutoff)
 Creates a transform plan.
void nfsft_precompute (int N, double kappa, unsigned int nfsft_flags, unsigned int fpt_flags)
 Performes precomputation up to the next power of two with respect to a given bandwidth $N \in \mathbb{N}_2$ .
void nfsft_forget ()
 Forgets all precomputed data.
void ndsft_trafo (nfsft_plan *plan)
 Executes a direct NDSFT, i.e.
void ndsft_adjoint (nfsft_plan *plan)
 Executes a direct adjoint NDSFT, i.e.
void nfsft_trafo (nfsft_plan *plan)
 Executes a NFSFT, i.e.
void nfsft_adjoint (nfsft_plan *plan)
 Executes an adjoint NFSFT, i.e.
void nfsft_finalize (nfsft_plan *plan)
 Destroys a plan.
void nfsft_precompute_x (nfsft_plan *plan)
fpt_set fpt_init (const int M, const int t, const unsigned int flags)
 Initializes a set of precomputed data for DPT transforms of equal length.
void fpt_precompute (fpt_set set, const int m, const double *alpha, const double *beta, const double *gamma, int k_start, const double threshold)
 Computes the data required for a single DPT transform.
void dpt_trafo (fpt_set set, const int m, const double complex *x, double complex *y, const int k_end, const unsigned int flags)
 Computes a single DPT transform.
void fpt_trafo (fpt_set set, const int m, const double complex *x, double complex *y, const int k_end, const unsigned int flags)
 Computes a single DPT transform.
void dpt_transposed (fpt_set set, const int m, double complex *x, double complex *y, const int k_end, const unsigned int flags)
 Computes a single DPT transform.
void fpt_transposed (fpt_set set, const int m, double complex *x, const double complex *y, const int k_end, const unsigned int flags)
 Computes a single DPT transform.
void fpt_finalize (fpt_set set)
void infft_init (infft_plan *ths, nfft_plan *mv)
 Simple initialisation.
void infft_init_advanced (infft_plan *ths, nfft_plan *mv, unsigned infft_flags)
 Advanced initialisation.
void infft_before_loop (infft_plan *ths)
 Setting up residuals before the actual iteration.
void infft_loop_one_step (infft_plan *ths)
 Doing one step in the iteration.
void infft_finalize (infft_plan *ths)
 Destroys the plan for the inverse transform.
void infct_init (infct_plan *ths, nfct_plan *mv)
 Simple initialisation.
void infct_init_advanced (infct_plan *ths, nfct_plan *mv, unsigned infct_flags)
 Advanced initialisation.
void infct_before_loop (infct_plan *ths)
 Setting up residuals before the actual iteration.
void infct_loop_one_step (infct_plan *ths)
 Doing one step in the iteration.
void infct_finalize (infct_plan *ths)
 Destroys the plan for the inverse transform.
void infst_init (infst_plan *ths, nfst_plan *mv)
 Simple initialisation.
void infst_init_advanced (infst_plan *ths, nfst_plan *mv, unsigned infst_flags)
 Advanced initialisation.
void infst_before_loop (infst_plan *ths)
 Setting up residuals before the actual iteration.
void infst_loop_one_step (infst_plan *ths)
 Doing one step in the iteration.
void infst_finalize (infst_plan *ths)
 Destroys the plan for the inverse transform.
void innfft_init (innfft_plan *ths, nnfft_plan *mv)
 Simple initialisation.
void innfft_init_advanced (innfft_plan *ths, nnfft_plan *mv, unsigned innfft_flags)
 Advanced initialisation.
void innfft_before_loop (innfft_plan *ths)
 Setting up residuals before the actual iteration.
void innfft_loop_one_step (innfft_plan *ths)
 Doing one step in the iteration.
void innfft_finalize (innfft_plan *ths)
 Destroys the plan for the inverse transform.
void imri_inh_2d1d_init (imri_inh_2d1d_plan *ths, mri_inh_2d1d_plan *mv)
 Simple initialisation.
void imri_inh_2d1d_init_advanced (imri_inh_2d1d_plan *ths, mri_inh_2d1d_plan *mv, unsigned imri_inh_2d1d_flags)
 Advanced initialisation.
void imri_inh_2d1d_before_loop (imri_inh_2d1d_plan *ths)
 Setting up residuals before the actual iteration.
void imri_inh_2d1d_loop_one_step (imri_inh_2d1d_plan *ths)
 Doing one step in the iteration.
void imri_inh_2d1d_finalize (imri_inh_2d1d_plan *ths)
 Destroys the plan for the inverse transform.
void imri_inh_3d_init (imri_inh_3d_plan *ths, mri_inh_3d_plan *mv)
 Simple initialisation.
void imri_inh_3d_init_advanced (imri_inh_3d_plan *ths, mri_inh_3d_plan *mv, unsigned imri_inh_3d_flags)
 Advanced initialisation.
void imri_inh_3d_before_loop (imri_inh_3d_plan *ths)
 Setting up residuals before the actual iteration.
void imri_inh_3d_loop_one_step (imri_inh_3d_plan *ths)
 Doing one step in the iteration.
void imri_inh_3d_finalize (imri_inh_3d_plan *ths)
 Destroys the plan for the inverse transform.
void infsft_init (infsft_plan *ths, nfsft_plan *mv)
 Simple initialisation.
void infsft_init_advanced (infsft_plan *ths, nfsft_plan *mv, unsigned infsft_flags)
 Advanced initialisation.
void infsft_before_loop (infsft_plan *ths)
 Setting up residuals before the actual iteration.
void infsft_loop_one_step (infsft_plan *ths)
 Doing one step in the iteration.
void infsft_finalize (infsft_plan *ths)
 Destroys the plan for the inverse transform.


Detailed Description

Header file for the nfft3 library.

Definition in file nfft3.h.


Define Documentation

#define MACRO_MV_PLAN float_type   ) 
 

Value:

int N_total;                          \
  int M_total;                          \
                                                                              \
  float_type *f_hat;                    \
  float_type *f;
Macros for public members inherited by all plan structures. Vector of samples, \ size is M_total float types.

Definition at line 14 of file nfft3.h.


Generated on 1 Nov 2006 by Doxygen 1.4.4