Main Page | Modules | Data Structures | Directories | File List | Data Fields | Globals | Related Pages

nfft3.h

Go to the documentation of this file.
00001 
00004 #ifndef NFFT3_H
00005 #define NFFT3_H
00006 
00008 #include <complex.h>
00010 #include <fftw3.h>
00011 
00013 #define MACRO_MV_PLAN(float_type)                                           \
00014 int N_total;                          \
00015 int M_total;                          \
00016 float_type *f_hat;                    \
00017 float_type *f;                        \
00018 
00019 /*###########################################################################*/
00020 /*###########################################################################*/
00021 /*###########################################################################*/
00022 
00032 #define PRE_PHI_HUT      (1U<< 0)
00033 #define FG_PSI           (1U<< 1)
00034 #define PRE_LIN_PSI      (1U<< 2)
00035 #define PRE_FG_PSI       (1U<< 3)
00036 #define PRE_PSI          (1U<< 4)
00037 #define PRE_FULL_PSI     (1U<< 5)
00038 #define MALLOC_X         (1U<< 6)
00039 #define MALLOC_F_HAT     (1U<< 7)
00040 #define MALLOC_F         (1U<< 8)
00041 #define FFT_OUT_OF_PLACE (1U<< 9)
00042 #define FFTW_INIT        (1U<< 10)
00043 
00044 #define MALLOC_V         (1U<< 11)
00045 
00046 #define SNDFT            (1U<< 12)
00047 
00048 #define PRE_ONE_PSI (PRE_LIN_PSI| PRE_FG_PSI| PRE_PSI| PRE_FULL_PSI)
00049 
00050 typedef struct nfft_plan_
00051 {
00053   MACRO_MV_PLAN(complex);
00054 
00055   int d;                                
00056   int *N;                               
00057   double *sigma;                        
00058   int *n;                               
00059   int n_total;                          
00060   int m;                                
00061   double *b;                            
00062   int K;                                
00064   unsigned nfft_flags;                  
00065   unsigned fftw_flags;                  
00067   double *x;                            
00069   double MEASURE_TIME_t[3];             
00072   fftw_plan  my_fftw_plan1;             
00073   fftw_plan  my_fftw_plan2;             
00075   double **c_phi_inv;                   
00076   double *psi;                          
00077   int *psi_index_g;                     
00078   int *psi_index_f;                     
00080   complex *g;
00081   complex *g_hat;
00082   complex *g1;                          
00083   complex *g2;                          
00085   double *spline_coeffs;            
00087 } nfft_plan;
00088 
00089 
00099 void ndft_trafo(nfft_plan *ths);
00100 
00110 void ndft_adjoint(nfft_plan *ths);
00111 
00122 void nfft_trafo(nfft_plan *ths);
00123 
00134 void nfft_adjoint(nfft_plan *ths);
00135 
00145 void nfft_init_1d(nfft_plan *ths, int N1, int M);
00146 
00157 void nfft_init_2d(nfft_plan *ths, int N1, int N2, int M);
00158 
00170 void nfft_init_3d(nfft_plan *ths, int N1, int N2, int N3, int M);
00171 
00182 void nfft_init(nfft_plan *ths, int d, int *N, int M);
00183 
00196 void nfft_init_advanced(nfft_plan *ths, int d, int *N, int M,
00197       unsigned nfft_flags_on, unsigned nfft_flags_off);
00198 
00213 void nfft_init_guru(nfft_plan *ths, int d, int *N, int M, int *n,
00214               int m, unsigned nfft_flags, unsigned fftw_flags);
00215 
00216 void nfft_precompute_full_psi(nfft_plan *ths);
00217 
00230 void nfft_precompute_one_psi(nfft_plan *ths);
00231 
00239 void nfft_finalize(nfft_plan *ths);
00240 
00244 /*###########################################################################*/
00245 /*###########################################################################*/
00246 /*###########################################################################*/
00247 
00255 typedef struct nfct_plan_
00256 {
00258   MACRO_MV_PLAN(double);
00259 
00260   int d;                                
00261   int *N;                               
00262   int *n;                               
00263   double *sigma;                        
00264   int m;                                
00266   double nfct_full_psi_eps;
00267   double *b;                            
00269   unsigned nfct_flags;                  
00270   unsigned fftw_flags;                  
00272   double *x;                            
00274   double MEASURE_TIME_t[3];             
00277   fftw_plan  my_fftw_r2r_plan;          
00278   fftw_r2r_kind *r2r_kind;              
00280   double **c_phi_inv;                   
00281   double *psi;                          
00282   int size_psi;                         
00283   int *psi_index_g;                     
00284   int *psi_index_f;                     
00286   double *g;
00287   double *g_hat;
00288   double *g1;                           
00289   double *g2;                           
00291   double *spline_coeffs;                
00293 } nfct_plan;
00294 
00295 
00305 void nfct_init_1d( nfct_plan *ths_plan, int N0, int M_total);
00306 
00317 void nfct_init_2d( nfct_plan *ths_plan, int N0, int N1, int M_total);
00318 
00330 void nfct_init_3d( nfct_plan *ths_plan, int N0, int N1, int N2, int M_total);
00331 
00342 void nfct_init( nfct_plan *ths_plan, int d, int *N, int M_total);
00343 
00358 void nfct_init_guru( nfct_plan *ths_plan, int d, int *N, int M_total, int *n,
00359                          int m, unsigned nfct_flags, unsigned fftw_flags);
00360 
00370 void nfct_precompute_psi( nfct_plan *ths_plan);
00371 
00372 
00381 void nfct_trafo( nfct_plan *ths_plan);
00382 
00391 void ndct_trafo( nfct_plan *ths_plan);
00392 
00401 void nfct_adjoint( nfct_plan *ths_plan);
00402 
00411 void ndct_adjoint( nfct_plan *ths_plan);
00412 
00420 void nfct_finalize( nfct_plan *ths_plan);
00421 
00431 double nfct_phi_hut( nfct_plan *ths_plan, int k, int d);
00432 
00442 double nfct_phi ( nfct_plan *ths_plan, double x, int d);
00443 
00451 int nfct_fftw_2N( int n);
00452 
00460 int nfct_fftw_2N_rev( int n);
00461 
00465 /*###########################################################################*/
00466 /*###########################################################################*/
00467 /*###########################################################################*/
00468 
00476 typedef struct nfst_plan_
00477 {
00479   MACRO_MV_PLAN(double);
00480 
00481   int d;                                
00482   int *N;                               
00483   int *n;                               
00484   double *sigma;                        
00485   int m;                                
00487   double nfst_full_psi_eps;
00488   double *b;                            
00490   unsigned nfst_flags;                  
00491   unsigned fftw_flags;                  
00493   double *x;                            
00495   double MEASURE_TIME_t[3];             
00498   fftw_plan  my_fftw_r2r_plan;         
00499   fftw_r2r_kind *r2r_kind;              
00501   double **c_phi_inv;                   
00502   double *psi;                          
00503   int size_psi;                         
00504   int *psi_index_g;                     
00505   int *psi_index_f;                     
00508   double *g;
00509   double *g_hat;
00510   double *g1;                           
00511   double *g2;                           
00513   double *spline_coeffs;                
00515 } nfst_plan;
00516 
00517 
00527 void nfst_init_1d( nfst_plan *ths_plan, int N0, int M_total);
00528 
00539 void nfst_init_2d( nfst_plan *ths_plan, int N0, int N1, int M_total);
00540 
00552 void nfst_init_3d( nfst_plan *ths_plan, int N0, int N1, int N2, int M_total);
00553 
00564 void nfst_init( nfst_plan *ths_plan, int d, int *N, int M_total);
00565 
00578 void nfst_init_m( nfst_plan *ths_plan, int d, int *N, int M_total, int m);
00579 
00594 void nfst_init_guru( nfst_plan *ths_plan, int d, int *N, int M_total, int *n,
00595                      int m, unsigned nfst_flags, unsigned fftw_flags);
00596 
00606 void nfst_precompute_psi( nfst_plan *ths_plan);
00607 
00616 void nfst_trafo( nfst_plan *ths_plan);
00617 
00626 void ndst_trafo( nfst_plan *ths_plan);
00627 
00628 
00629 
00638 void nfst_adjoint( nfst_plan *ths_plan);
00639 
00648 void ndst_adjoint( nfst_plan *ths_plan);
00649 
00657 void nfst_finalize( nfst_plan *ths_plan);
00658 
00665 void nfst_full_psi( nfst_plan *ths_plan, double eps);
00666 
00676 double nfst_phi_hut( nfst_plan *ths_plan, int k, int d);
00677 
00687 double nfst_phi ( nfst_plan *ths_plan, double x, int d);
00688 
00696 int nfst_fftw_2N( int n);
00697 
00705 int nfst_fftw_2N_rev( int n);
00706 
00710 /*###########################################################################*/
00711 /*###########################################################################*/
00712 /*###########################################################################*/
00713 
00723 typedef struct
00724 {
00726   MACRO_MV_PLAN(complex);
00727 
00728   int d;                                
00729   double *sigma;                        
00730   double *a;                            
00731   int *N;                               
00732   int *N1;                              
00733   int *aN1;                             
00734   int m;                                
00735   double *b;                            
00736   int K;                                
00738   int aN1_total;                        
00740   nfft_plan *direct_plan;               
00741   unsigned nnfft_flags;                 
00742   int *n;                               
00744   double *x;                            
00745   double *v;                            
00747   double *c_phi_inv;                    
00748   double *psi;                          
00749   int size_psi;                         
00750   int *psi_index_g;                     
00751   int *psi_index_f;                     
00752   complex *F;
00753 
00754   double *spline_coeffs;                
00756 } nnfft_plan;
00757 
00758 
00770 void nnfft_init(nnfft_plan *ths_plan, int d, int N_total, int M_total, int *N);
00771 
00786 void nnfft_init_guru(nnfft_plan *ths_plan, int d, int N_total, int M_total,
00787                      int *N, int *N1, int m, unsigned nnfft_flags);
00788 
00799 void nndft_trafo(nnfft_plan *ths_plan);
00800 
00811 void nndft_adjoint(nnfft_plan *ths_plan);
00812 
00823 void nnfft_trafo(nnfft_plan *ths_plan);
00824 
00835 void nnfft_adjoint(nnfft_plan *ths_plan);
00836 
00848 void nnfft_precompute_lin_psi(nnfft_plan *ths_plan);
00849 
00862 void nnfft_precompute_psi(nnfft_plan *ths_plan);
00863 
00877 void nnfft_precompute_full_psi(nnfft_plan *ths_plan);
00878 
00891 void nnfft_precompute_phi_hut(nnfft_plan *ths_plan);
00892 
00900 void nnfft_finalize(nnfft_plan *ths_plan);
00901 
00905 /*###########################################################################*/
00906 /*###########################################################################*/
00907 /*###########################################################################*/
00908 
00914 typedef struct nsfft_plan_
00915 {
00916   MACRO_MV_PLAN(complex);
00917 
00918   int d;                                
00919   int J;                                
00923   int sigma;                            
00925   unsigned flags;                       
00927   int *index_sparse_to_full;            
00930   int r_act_nfft_plan;                  
00931   nfft_plan *act_nfft_plan;             
00932   nfft_plan *center_nfft_plan;          
00934   fftw_plan* set_fftw_plan1;            
00935   fftw_plan* set_fftw_plan2;            
00937   nfft_plan *set_nfft_plan_1d;          
00938   nfft_plan *set_nfft_plan_2d;          
00940   double *x_transposed;                 
00941   double *x_102,*x_201,*x_120,*x_021;   
00943 } nsfft_plan;
00944 
00954 void nsdft_trafo(nsfft_plan *ths);
00955 
00965 void nsdft_adjoint(nsfft_plan *ths);
00966 
00976 void nsfft_trafo(nsfft_plan *ths);
00977 
00987 void nsfft_adjoint(nsfft_plan *ths);
00988 
00996 void nsfft_cp(nsfft_plan *ths, nfft_plan *ths_nfft);
00997 
01005 void nsfft_init_random_nodes_coeffs(nsfft_plan *ths);
01006 
01019 void nsfft_init(nsfft_plan *ths, int d, int J, int M, int m, unsigned flags);
01020 
01028 void nsfft_finalize(nsfft_plan *ths);
01029 
01033 /*###########################################################################*/
01034 /*###########################################################################*/
01035 /*###########################################################################*/
01036 
01044 typedef struct
01045 {
01047   MACRO_MV_PLAN(complex);
01048 
01049   nfft_plan plan;
01050 
01051   int N3;
01052   double sigma3;
01053   double *t;
01054   double *w;
01055 } mri_inh_2d1d_plan;
01056 
01060 typedef struct
01061 {
01063   MACRO_MV_PLAN(complex);
01064 
01065   nfft_plan plan;
01066 
01067   int N3;
01068   double sigma3;
01069   double *t;
01070   double *w;
01071 } mri_inh_3d_plan;
01072 
01073 
01086 void mri_inh_2d1d_trafo(mri_inh_2d1d_plan *ths);
01087 
01100 void mri_inh_2d1d_adjoint(mri_inh_2d1d_plan *ths);
01101 
01115 void mri_inh_2d1d_init_guru(mri_inh_2d1d_plan *ths, int *N, int M, int *n,
01116                     int m, double sigma, unsigned nfft_flags, unsigned fftw_flags);
01117 
01125 void mri_inh_2d1d_finalize(mri_inh_2d1d_plan *ths);
01126 
01139 void mri_inh_3d_trafo(mri_inh_3d_plan *ths);
01140 
01153 void mri_inh_3d_adjoint(mri_inh_3d_plan *ths);
01154 
01155 void mri_inh_3d_init_guru(mri_inh_3d_plan *ths, int *N, int M, int *n,
01156                     int m, double sigma, unsigned nfft_flags, unsigned fftw_flags);
01157 
01165 void mri_inh_3d_finalize(mri_inh_3d_plan *ths);
01166 
01170 /*###########################################################################*/
01171 /*###########################################################################*/
01172 /*###########################################################################*/
01173 
01347 #define TEXTURE_MAX_ANGLE (2*3.1415926535897932384)
01348 
01357 #define TEXTURE_DEF_PRECOMPUTE_FLAGS 0U
01358 
01362 #define TEXTURE_DEF_NFSFT_PRECOMPUTE_FLAGS 0U
01363 
01367 #define TEXTURE_DEF_NFSFT_THRESHOLD 1000.0
01368 
01372 #define TEXTURE_DEF_INIT_FLAGS 0U
01373 
01377 #define TEXTURE_DEF_NFSFT_INIT_FLAGS 0U
01378 
01382 #define TEXTURE_DEF_NFFT_CUTOFF 8
01383 
01400 typedef struct texture_plan_ {
01401 
01419   MACRO_MV_PLAN(complex);
01420 
01424   int N;
01425 
01429   int N1;
01430 
01434   int N2;
01435 
01439   const double *h_phi;
01440 
01444   const double *h_theta;
01445 
01449   const double *r;
01450 
01454   unsigned int nfsft_init_flags;
01455 
01459   unsigned int nfft_cutoff;
01460 
01463   double *cos_h_theta;
01464 
01467   double *sin_h_theta;
01468 
01471   complex **nfsft_f_hat;
01472 
01475   complex *nfsft_f;
01476 
01479   double *nfsft_angles;
01480 } texture_plan;
01481 
01494 void texture_precompute(int N);
01495 
01509 void texture_precompute_advanced(int N, unsigned int texture_precompute_flags,
01510     unsigned int nfsft_precompute_flags, double nfsft_threshold);
01511 
01540 void texture_init(texture_plan *ths, int N, int N1, int N2, complex* omega,
01541     complex* x, const double* h_phi, const double* h_theta, const double* r);
01542 
01578 void texture_init_advanced(texture_plan *ths, int N, int N1, int N2,
01579     complex* omega, complex* x, const double* h_phi, const double* h_theta,
01580     const double *r, unsigned int texture_init_flags,
01581     unsigned int nfsft_init_flags, int nfft_cutoff);
01582 
01595 void texture_trafo(texture_plan *ths);
01596 
01609 void texture_adjoint(texture_plan *ths);
01610 
01613 void texture_finalize(texture_plan *ths);
01614 
01618 void texture_forget();
01619 
01638 int texture_flat_index(int l, int m, int n);
01639 
01647 int texture_flat_length(int N);
01648 
01653 int texture_get_omega_length(texture_plan *ths);
01654 
01659 int texture_get_x_length(texture_plan *ths);
01660 
01665 int texture_get_N(texture_plan *ths);
01666 
01671 int texture_get_N1(texture_plan *ths);
01672 
01677 int texture_get_N2(texture_plan *ths);
01678 
01683 const complex *texture_get_omega(texture_plan *ths);
01684 
01691 void texture_set_omega(texture_plan *ths, complex* omega);
01692 
01697 const complex *texture_get_x(texture_plan *ths);
01698 
01705 void texture_set_x(texture_plan *ths, complex* x);
01706 
01711 const double *texture_get_h_phi(texture_plan *ths);
01712 
01719 void texture_set_h_phi(texture_plan *ths, const double* h_phi);
01720 
01725 const double *texture_get_h_theta(texture_plan *ths);
01726 
01733 void texture_set_h_theta(texture_plan *ths, const double* h_theta);
01734 
01739 const double *texture_get_r(texture_plan *ths);
01740 
01747 void texture_set_r(texture_plan *ths, const double* r);
01748 
01754 unsigned int texture_get_nfsft_init_flags(texture_plan *ths);
01755 
01761 void texture_set_nfsft_init_flags(texture_plan *ths,
01762     unsigned int nfsft_init_flags);
01763 
01769 int texture_get_nfft_cutoff(texture_plan *ths);
01770 
01777 void texture_set_nfft_cutoff(texture_plan *ths, int nfft_cutoff);
01778 
01785 /*###########################################################################*/
01786 /*###########################################################################*/
01787 /*###########################################################################*/
01788 
02040 /* Planner flags */
02041 
02061 #define NFSFT_NORMALIZED    (1U << 0)
02062 
02073 #define NFSFT_USE_NDFT      (1U << 1)
02074 
02086 #define NFSFT_USE_DPT       (1U << 2)
02087 
02101 #define NFSFT_MALLOC_X      (1U << 3)
02102 
02116 #define NFSFT_MALLOC_F_HAT  (1U << 5)
02117 
02131 #define NFSFT_MALLOC_F      (1U << 6)
02132 
02143 #define NFSFT_PRESERVE_F_HAT (1U << 7)
02144 
02156 #define NFSFT_PRESERVE_X     (1U << 8)
02157 
02168 #define NFSFT_PRESERVE_F     (1U << 9)
02169 
02179 #define NFSFT_DESTROY_F_HAT    (1U << 10)
02180 
02191 #define NFSFT_DESTROY_X      (1U << 11)
02192 
02202 #define NFSFT_DESTROY_F      (1U << 12)
02203 
02204 /* Precomputation flags */
02205 
02215 #define NFSFT_NO_DIRECT_ALGORITHM    (1U << 13)
02216 
02226 #define NFSFT_NO_FAST_ALGORITHM      (1U << 14)
02227 
02235 #define NFSFT_ZERO_F_HAT             (1U << 16)
02236 
02237 /* */
02238 
02247 #define NFSFT_INDEX(k,n,plan)        ((2*(plan)->N+2)*((plan)->N-n+1)+(plan)->N+k+1)
02248 
02253 #define NFSFT_F_HAT_SIZE(N)          ((2*N+2)*(2*N+2))
02254 
02256 typedef struct nfsft_plan_
02257 {
02259   MACRO_MV_PLAN(complex);
02260 
02261   /* Public members */
02262   int N;                              
02263   double *x;                          
02270   /* Private members */
02271   /*int NPT;*/                        
02273   int t;                              
02275   unsigned int flags;                 
02276   nfft_plan plan_nfft;                
02277   complex *f_hat_intern;              
02279 } nfsft_plan;
02280 
02290 void nfsft_init(nfsft_plan *plan, int N, int M);
02291 
02302 void nfsft_init_advanced(nfsft_plan* plan, int N, int M, unsigned int
02303                          nfsft_flags);
02304 
02316 void nfsft_init_guru(nfsft_plan *plan, int N, int M, unsigned int nfsft_flags,
02317   int nfft_flags, int nfft_cutoff);
02318 
02331 void nfsft_precompute(int N, double kappa, unsigned int nfsft_flags,
02332   unsigned int fpt_flags);
02333 
02339 void nfsft_forget();
02340 
02352 void ndsft_trafo(nfsft_plan* plan);
02353 
02366 void ndsft_adjoint(nfsft_plan* plan);
02367 
02379 void nfsft_trafo(nfsft_plan* plan);
02380 
02393 void nfsft_adjoint(nfsft_plan* plan);
02394 
02402 void nfsft_finalize(nfsft_plan* plan);
02403 
02404 void nfsft_precompute_x(nfsft_plan *plan);
02405 
02406 
02410 /*###########################################################################*/
02411 /*###########################################################################*/
02412 /*###########################################################################*/
02413 
02422 /* Flags for fpt_init() */
02423 #define FPT_NO_FAST_ALGORITHM (1U << 2) 
02424 #define FPT_NO_DIRECT_ALGORITHM (1U << 3) 
02425 #define FPT_NO_STABILIZATION  (1U << 0) 
02428 #define FPT_PERSISTENT_DATA   (1U << 4) 
02430 /* Flags for fpt_trafo(), dpt_transposed(), fpt_trafo(), fpt_transposed() */
02431 #define FPT_FUNCTION_VALUES   (1U << 5) 
02434 #define FPT_AL_SYMMETRY       (1U << 6) 
02436 /* Data structures */
02437 typedef struct fpt_set_s_ *fpt_set;     
02454 fpt_set fpt_init(const int M, const int t, const unsigned int flags);
02455 
02476 void fpt_precompute(fpt_set set, const int m, const double *alpha,
02477                     const double *beta, const double *gamma, int k_start,
02478                     const double threshold);
02479 
02490 void dpt_trafo(fpt_set set, const int m, const complex *x, complex *y,
02491   const int k_end, const unsigned int flags);
02492 
02503 void fpt_trafo(fpt_set set, const int m, const complex *x, complex *y,
02504   const int k_end, const unsigned int flags);
02505 
02516 void dpt_transposed(fpt_set set, const int m, complex *x, const complex *y,
02517   const int k_end, const unsigned int flags);
02518 
02529 void fpt_transposed(fpt_set set, const int m, complex *x, const complex *y,
02530   const int k_end, const unsigned int flags);
02531 
02532 void fpt_finalize(fpt_set set);
02533 
02537 /*###########################################################################*/
02538 /*###########################################################################*/
02539 /*###########################################################################*/
02540 
02548 #define LANDWEBER             (1U<< 0)
02549 #define STEEPEST_DESCENT      (1U<< 1)
02550 #define CGNR                  (1U<< 2)
02551 #define CGNE                  (1U<< 3)
02552 #define NORMS_FOR_LANDWEBER   (1U<< 4)
02553 #define PRECOMPUTE_WEIGHT     (1U<< 5)
02554 #define PRECOMPUTE_DAMP       (1U<< 6)
02555 
02571 #define F(MV, FLT, name, ...) void i ## MV ## _ ## name(__VA_ARGS__)
02572 
02573 #define MACRO_SOLVER_PLAN(MV, FLT)                                            \
02574 typedef struct i ## MV ## _plan_                                              \
02575 {                                                                        \
02576   MV ## _plan *mv;                      \
02577   unsigned flags;      \
02578                                                                               \
02579   double *w;                           \
02580   double *w_hat;                       \
02581                                                                               \
02582   FLT *y;                               \
02583                                                                               \
02584   FLT *f_hat_iter;                      \
02585                         \
02586   FLT *r_iter;              \
02587   FLT *z_hat_iter;                       \
02588   FLT *p_hat_iter;                       \
02589   FLT *v_iter;                           \
02590                         \
02591   double alpha_iter;                    \
02592   double beta_iter;                     \
02593                         \
02594   double dot_r_iter;                    \
02595   double dot_r_iter_old;                \
02596   double dot_z_hat_iter;                \
02597   double dot_z_hat_iter_old;            \
02598   double dot_p_hat_iter;                \
02599   double dot_v_iter;                    \
02600 } i ## MV ## _plan;                                                           \
02601                         \
02602 F(MV, FLT, init,    i ## MV ## _plan *ths, MV ## _plan *mv);            \
02603 F(MV, FLT, init_advanced, i ## MV ## _plan *ths, MV ## _plan *mv,             \
02604              unsigned i ## MV ## _flags);                  \
02605 F(MV, FLT, before_loop,   i ## MV ## _plan *ths);                             \
02606 F(MV, FLT, loop_one_step, i ## MV ## _plan *ths);                             \
02607 F(MV, FLT, finalize,      i ## MV ## _plan *ths);                             \
02608 
02609 
02610 MACRO_SOLVER_PLAN(nfft, complex)
02611 MACRO_SOLVER_PLAN(nfct, double)
02612 MACRO_SOLVER_PLAN(nfst, double)
02613 MACRO_SOLVER_PLAN(nnfft, complex)
02614 MACRO_SOLVER_PLAN(mri_inh_2d1d, complex)
02615 MACRO_SOLVER_PLAN(mri_inh_3d, complex)
02616 MACRO_SOLVER_PLAN(nfsft, complex)
02617 MACRO_SOLVER_PLAN(texture, complex)
02621 #endif
02622 /* nfft3.h */

Generated on Wed May 10 20:39:39 2006 for NFFT by  doxygen 1.4.4