NFFT Logo 3.1.1 API Reference

nnfft.c

00001 /*
00002  * Copyright (c) 2002, 2009 Jens Keiner, Stefan Kunis, Daniel Potts
00003  *
00004  * This program is free software; you can redistribute it and/or modify it under
00005  * the terms of the GNU General Public License as published by the Free Software
00006  * Foundation; either version 2 of the License, or (at your option) any later
00007  * version.
00008  *
00009  * This program is distributed in the hope that it will be useful, but WITHOUT
00010  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00011  * FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
00012  * details.
00013  *
00014  * You should have received a copy of the GNU General Public License along with
00015  * this program; if not, write to the Free Software Foundation, Inc., 51
00016  * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
00017  */
00018 
00019 /* $Id: nnfft.c 3198 2009-05-27 14:16:50Z keiner $ */
00020 
00021 #include <stdio.h>
00022 #include <math.h>
00023 #include <string.h>
00024 #include <stdlib.h>
00025 
00026 #include <complex.h>
00027 
00028 #include "nfft3util.h"
00029 #include "nfft3.h"
00030 
00031 
00032 #define MACRO_nndft_init_result_trafo memset(f,0,ths->M_total*sizeof(double _Complex));
00033 #define MACRO_nndft_init_result_conjugated MACRO_nndft_init_result_trafo
00034 #define MACRO_nndft_init_result_adjoint memset(f_hat,0,ths->N_total*sizeof(double _Complex));
00035 #define MACRO_nndft_init_result_transposed MACRO_nndft_init_result_adjoint
00036 
00037 #define MACRO_nndft_sign_trafo      (-2.0*PI)
00038 #define MACRO_nndft_sign_conjugated (+2.0*PI)
00039 #define MACRO_nndft_sign_adjoint    (+2.0*PI)
00040 #define MACRO_nndft_sign_transposed (-2.0*PI)
00041 
00042 #define MACRO_nndft_compute_trafo (*fj) += (*f_hat_k)*cexp(+ _Complex_I*omega);
00043 
00044 #define MACRO_nndft_compute_conjugated MACRO_nndft_compute_trafo
00045 
00046 #define MACRO_nndft_compute_adjoint (*f_hat_k) += (*fj)*cexp(+ _Complex_I*omega);
00047 
00048 #define MACRO_nndft_compute_transposed MACRO_nndft_compute_adjoint
00049 
00050 #define MACRO_nndft(which_one)                                                \
00051 void nndft_ ## which_one (nnfft_plan *ths)                                    \
00052 {                                                                             \
00053   int j;                               \
00054   int t;                               \
00055   int l;                               \
00056   double _Complex *f_hat, *f;          \
00057   double _Complex *f_hat_k;            \
00058   double _Complex *fj;                 \
00059   double omega;                        \
00060                                                                               \
00061   f_hat=ths->f_hat; f=ths->f;                                                 \
00062                                                                               \
00063   MACRO_nndft_init_result_ ## which_one                                       \
00064                                                                               \
00065   for(j=0, fj=f; j<ths->M_total; j++, fj++)                                   \
00066   {                                                                           \
00067     for(l=0, f_hat_k=f_hat; l<ths->N_total; l++, f_hat_k++)                   \
00068     {                                                                         \
00069       omega=0.0;                                                              \
00070       for(t = 0; t<ths->d; t++)                                               \
00071         omega+=ths->v[l*ths->d+t] * ths->x[j*ths->d+t] * ths->N[t];           \
00072                                                                               \
00073       omega*= MACRO_nndft_sign_ ## which_one;                                 \
00074                                                                               \
00075       MACRO_nndft_compute_ ## which_one                                       \
00076                                                                               \
00077      } /* for(l) */                                                           \
00078    } /* for(j) */                                                             \
00079 } /* nndft_trafo */                                                           \
00080 
00081 MACRO_nndft(trafo)
00082 MACRO_nndft(adjoint)
00083 
00086 void nnfft_uo(nnfft_plan *ths,int j,int *up,int *op,int act_dim)
00087 {
00088   double c;
00089   int u,o;
00090 
00091   c = ths->v[j*ths->d+act_dim] * ths->n[act_dim];
00092 
00093   u = c; o = c;
00094   if(c < 0)
00095     u = u-1;
00096   else
00097     o = o+1;
00098 
00099   u = u - (ths->m); o = o + (ths->m);
00100 
00101   up[0]=u; op[0]=o;
00102 }
00103 
00107 #define MACRO_nnfft_B_init_result_A memset(f,0,ths->N_total*sizeof(double _Complex));
00108 #define MACRO_nnfft_B_init_result_T memset(g,0,ths->aN1_total*sizeof(double _Complex));
00109 
00110 #define MACRO_nnfft_B_PRE_FULL_PSI_compute_A {                                \
00111   (*fj) += ths->psi[ix] * g[ths->psi_index_g[ix]];                            \
00112 }
00113 
00114 #define MACRO_nnfft_B_PRE_FULL_PSI_compute_T {                                \
00115   g[ths->psi_index_g[ix]] += ths->psi[ix] * (*fj);                            \
00116 }
00117 
00118 #define MACRO_nnfft_B_compute_A {                                             \
00119   (*fj) += phi_prod[ths->d] * g[ll_plain[ths->d]];                            \
00120 }
00121 
00122 #define MACRO_nnfft_B_compute_T {                                             \
00123   g[ll_plain[ths->d]] += phi_prod[ths->d] * (*fj);                            \
00124 }
00125 
00126 #define MACRO_with_PRE_LIN_PSI (ths->psi[(ths->K+1)*t2+y_u[t2]]*              \
00127                                 (y_u[t2]+1-y[t2]) +                           \
00128                                 ths->psi[(ths->K+1)*t2+y_u[t2]+1]*            \
00129                                 (y[t2]-y_u[t2]))
00130 #define MACRO_with_PRE_PSI     ths->psi[(j*ths->d+t2)*(2*ths->m+2)+lj[t2]]
00131 #define MACRO_without_PRE_PSI  PHI(-ths->v[j*ths->d+t2]+                      \
00132                                ((double)l[t2])/ths->N1[t2], t2)
00133 
00134 #define MACRO_init_uo_l_lj_t {                                                \
00135   for(t = ths->d-1; t>=0; t--)                                                \
00136     {                                                                         \
00137       nnfft_uo(ths,j,&u[t],&o[t],t);                                          \
00138       l[t] = u[t];                                                            \
00139       lj[t] = 0;                                                              \
00140     } /* for(t) */                                                            \
00141   t++;                                                                        \
00142 }
00143 
00144 #define MACRO_update_with_PRE_PSI_LIN {                                       \
00145   for(t2=t; t2<ths->d; t2++)                                                  \
00146     {                                                                         \
00147       y[t2] = fabs(((-ths->N1[t2]*ths->v[j*ths->d+t2]+(double)l[t2])          \
00148           * ((double)ths->K))/(ths->m+1));                                    \
00149       y_u[t2] = (int)y[t2];                                                   \
00150     } /* for(t2) */                                                           \
00151 }
00152 
00153 #define MACRO_update_phi_prod_ll_plain(which_one) {                           \
00154   for(t2=t; t2<ths->d; t2++)                                                  \
00155     {                                                                         \
00156       phi_prod[t2+1]=phi_prod[t2]* MACRO_ ## which_one;                       \
00157       ll_plain[t2+1]=ll_plain[t2]*ths->aN1[t2] +                              \
00158                      (l[t2]+ths->aN1[t2]*3/2)%ths->aN1[t2];                   \
00159       /* 3/2 because of the (not needed) fftshift and to be in [0 aN1[t2]]?!*/\
00160     } /* for(t2) */                                                           \
00161 }
00162 
00163 #define MACRO_count_uo_l_lj_t {                                               \
00164   for(t = ths->d-1; (t>0)&&(l[t]==o[t]); t--)                                 \
00165     {                                                                         \
00166       l[t] = u[t];                                                            \
00167       lj[t] = 0;                                                              \
00168     } /* for(t) */                                                            \
00169                                                                               \
00170   l[t]++;                                                                     \
00171   lj[t]++;                                                                    \
00172 }
00173 
00174 #define MACRO_nnfft_B(which_one)                                              \
00175 static inline void nnfft_B_ ## which_one (nnfft_plan *ths)                    \
00176 {                                                                             \
00177   int lprod;                           \
00178   int u[ths->d], o[ths->d];            \
00179   int t, t2;                           \
00180   int j;                               \
00181   int l_L, ix;                         \
00182   int l[ths->d];                       \
00183   int lj[ths->d];                      \
00184   int ll_plain[ths->d+1];              \
00185   double phi_prod[ths->d+1];           \
00186   double _Complex *f, *g;              \
00187   double _Complex *fj;                 \
00188   double y[ths->d];                                                           \
00189   int y_u[ths->d];                                                            \
00190                                                                               \
00191   f=ths->f_hat; g=ths->F;                                                     \
00192                                                                               \
00193   MACRO_nnfft_B_init_result_ ## which_one                                     \
00194                                                                               \
00195   if(ths->nnfft_flags & PRE_FULL_PSI)                                         \
00196     {                                                                         \
00197       for(ix=0, j=0, fj=f; j<ths->N_total; j++,fj++)                          \
00198         for(l_L=0; l_L<ths->psi_index_f[j]; l_L++, ix++)                      \
00199           MACRO_nnfft_B_PRE_FULL_PSI_compute_ ## which_one;                   \
00200       return;                                                                 \
00201     }                                                                         \
00202                                                                               \
00203   phi_prod[0]=1;                                                              \
00204   ll_plain[0]=0;                                                              \
00205                                                                               \
00206   for(t=0,lprod = 1; t<ths->d; t++)                                           \
00207     lprod *= (2*ths->m+2);                                                    \
00208                                                                               \
00209   if(ths->nnfft_flags & PRE_PSI)                                              \
00210     {                                                                         \
00211       for(j=0, fj=f; j<ths->N_total; j++, fj++)                               \
00212         {                                                                     \
00213           MACRO_init_uo_l_lj_t;                                               \
00214                                                                               \
00215           for(l_L=0; l_L<lprod; l_L++)                                        \
00216             {                                                                 \
00217               MACRO_update_phi_prod_ll_plain(with_PRE_PSI);                   \
00218                                                                               \
00219               MACRO_nnfft_B_compute_ ## which_one;                            \
00220                                                                               \
00221               MACRO_count_uo_l_lj_t;                                          \
00222             } /* for(l_L) */                                                  \
00223         } /* for(j) */                                                        \
00224       return;                                                                 \
00225     } /* if(PRE_PSI) */                                                       \
00226                                                                               \
00227   if(ths->nnfft_flags & PRE_LIN_PSI)                                          \
00228     {                                                                         \
00229       for(j=0, fj=f; j<ths->N_total; j++, fj++)                               \
00230         {                                                                     \
00231           MACRO_init_uo_l_lj_t;                                               \
00232                                                                               \
00233           for(l_L=0; l_L<lprod; l_L++)                                        \
00234             {                                                                 \
00235               MACRO_update_with_PRE_PSI_LIN;                                  \
00236                                                                               \
00237               MACRO_update_phi_prod_ll_plain(with_PRE_LIN_PSI);               \
00238                                                                               \
00239               MACRO_nnfft_B_compute_ ## which_one;                            \
00240                                                                               \
00241               MACRO_count_uo_l_lj_t;                                          \
00242             } /* for(l_L) */                                                  \
00243         } /* for(j) */                                                        \
00244       return;                                                                 \
00245     } /* if(PRE_LIN_PSI) */                                                   \
00246                                                                               \
00247   /* no precomputed psi at all */                                             \
00248   for(j=0, fj=f; j<ths->N_total; j++, fj++)                                   \
00249     {                                                                         \
00250                                                                               \
00251       MACRO_init_uo_l_lj_t;                                                   \
00252                                                                               \
00253       for(l_L=0; l_L<lprod; l_L++)                                            \
00254         {                                                                     \
00255           MACRO_update_phi_prod_ll_plain(without_PRE_PSI);                    \
00256                                                                               \
00257           MACRO_nnfft_B_compute_ ## which_one;                                \
00258                                                                               \
00259           MACRO_count_uo_l_lj_t;                                              \
00260         } /* for(l_L) */                                                      \
00261     } /* for(j) */                                                            \
00262 } /* nnfft_B */
00263 
00264 MACRO_nnfft_B(A)
00265 MACRO_nnfft_B(T)
00266 
00267 static inline void nnfft_D (nnfft_plan *ths){
00268   int j,t;
00269   double tmp;
00270 
00271   if(ths->nnfft_flags & PRE_PHI_HUT)
00272   {
00273       for(j=0; j<ths->M_total; j++)
00274     ths->f[j] *= ths->c_phi_inv[j];
00275   }
00276   else
00277   {
00278       for(j=0; j<ths->M_total; j++)
00279       {
00280     tmp = 1.0;
00281     /* multiply with N1, because x was modified */
00282     for(t=0; t<ths->d; t++)
00283         tmp*= 1.0 /((PHI_HUT(ths->x[ths->d*j + t]*((double)ths->N[t]),t)) );
00284     ths->f[j] *= tmp;
00285       }
00286   }
00287 }
00288 
00291 void nnfft_trafo(nnfft_plan *ths)
00292 {
00293   int j,t;
00294 
00295   nnfft_B_T(ths);
00296 
00297   for(j=0;j<ths->M_total;j++) {
00298     for(t=0;t<ths->d;t++) {
00299       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
00300     }
00301   }
00302 
00303   /* allows for external swaps of ths->f */
00304   ths->direct_plan->f = ths->f;
00305 
00306   nfft_trafo(ths->direct_plan);
00307 
00308   for(j=0;j<ths->M_total;j++) {
00309     for(t=0;t<ths->d;t++) {
00310       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
00311     }
00312   }
00313 
00314   nnfft_D(ths);
00315 } /* nnfft_trafo */
00316 
00317 void nnfft_adjoint(nnfft_plan *ths)
00318 {
00319   int j,t;
00320 
00321   nnfft_D(ths);
00322 
00323   for(j=0;j<ths->M_total;j++) {
00324     for(t=0;t<ths->d;t++) {
00325       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
00326     }
00327   }
00328 
00329   /* allows for external swaps of ths->f */
00330   ths->direct_plan->f=ths->f;
00331 
00332   nfft_adjoint(ths->direct_plan);
00333 
00334   for(j=0;j<ths->M_total;j++) {
00335     for(t=0;t<ths->d;t++) {
00336       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
00337     }
00338   }
00339 
00340   nnfft_B_A(ths);
00341 } /* nnfft_adjoint */
00342 
00345 void nnfft_precompute_phi_hut(nnfft_plan *ths)
00346 {
00347   int j;                                
00348   int t;                                
00349   double tmp;
00350 
00351   ths->c_phi_inv= (double*)nfft_malloc(ths->M_total*sizeof(double));
00352 
00353   for(j=0; j<ths->M_total; j++)
00354     {
00355       tmp = 1.0;
00356       for(t=0; t<ths->d; t++)
00357         tmp*= 1.0 /(PHI_HUT(ths->x[ths->d*j + t]*((double)ths->N[t]),t));
00358       ths->c_phi_inv[j]=tmp;
00359     }
00360 } /* nnfft_phi_hut */
00361 
00362 
00365 void nnfft_precompute_lin_psi(nnfft_plan *ths)
00366 {
00367   int t;                                
00368   int j;                                
00369   double step;                          
00371   nfft_precompute_lin_psi(ths->direct_plan);
00372 
00373   for (t=0; t<ths->d; t++)
00374     {
00375       step=((double)(ths->m+1))/(ths->K*ths->N1[t]);
00376       for(j=0;j<=ths->K;j++)
00377         {
00378           ths->psi[(ths->K+1)*t + j] = PHI(j*step,t);
00379         } /* for(j) */
00380     } /* for(t) */
00381 }
00382 
00383 void nnfft_precompute_psi(nnfft_plan *ths)
00384 {
00385   int t;                                
00386   int j;                                
00387   int l;                                
00388   int lj;                               
00389   int u, o;                             
00391   for (t=0; t<ths->d; t++)
00392     for(j=0;j<ths->N_total;j++)
00393       {
00394         nnfft_uo(ths,j,&u,&o,t);
00395 
00396         for(l=u, lj=0; l <= o; l++, lj++)
00397           ths->psi[(j*ths->d+t)*(2*ths->m+2)+lj]=
00398             (PHI((-ths->v[j*ths->d+t]+((double)l)/((double)ths->N1[t])),t));
00399       } /* for(j) */
00400 
00401   for(j=0;j<ths->M_total;j++) {
00402     for(t=0;t<ths->d;t++) {
00403       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
00404     }
00405   }
00406 
00407   nfft_precompute_psi(ths->direct_plan);
00408 
00409   for(j=0;j<ths->M_total;j++) {
00410     for(t=0;t<ths->d;t++) {
00411       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
00412     }
00413   }
00414   /* for(t) */
00415 } /* nfft_precompute_psi */
00416 
00417 
00418 
00422 void nnfft_precompute_full_psi(nnfft_plan *ths)
00423 {
00424   int t,t2;                             
00425   int j;                                
00426   int l_L;                              
00427   int l[ths->d];                       
00428   int lj[ths->d];                      
00429   int ll_plain[ths->d+1];              
00430   int lprod;                            
00431   int u[ths->d], o[ths->d];           
00433   double phi_prod[ths->d+1];
00434 
00435   int ix,ix_old;
00436 
00437   for(j=0;j<ths->M_total;j++) {
00438     for(t=0;t<ths->d;t++) {
00439       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
00440     }
00441   }
00442 
00443   nnfft_precompute_psi(ths);
00444 
00445   nfft_precompute_full_psi(ths->direct_plan);
00446 
00447   for(j=0;j<ths->M_total;j++) {
00448     for(t=0;t<ths->d;t++) {
00449       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
00450     }
00451   }
00452 
00453   phi_prod[0]=1;
00454   ll_plain[0]=0;
00455 
00456   for(t=0,lprod = 1; t<ths->d; t++)
00457     lprod *= 2*ths->m+2;
00458 
00459   for(j=0,ix=0,ix_old=0; j<ths->N_total; j++)
00460     {
00461       MACRO_init_uo_l_lj_t;
00462 
00463       for(l_L=0; l_L<lprod; l_L++, ix++)
00464         {
00465           MACRO_update_phi_prod_ll_plain(without_PRE_PSI);
00466 
00467           ths->psi_index_g[ix]=ll_plain[ths->d];
00468           ths->psi[ix]=phi_prod[ths->d];
00469 
00470           MACRO_count_uo_l_lj_t;
00471         } /* for(l_L) */
00472 
00473 
00474       ths->psi_index_f[j]=ix-ix_old;
00475       ix_old=ix;
00476     } /* for(j) */
00477 }
00478 
00479 void nnfft_init_help(nnfft_plan *ths, int m2, unsigned nfft_flags, unsigned fftw_flags)
00480 {
00481   int t;                                
00482   int lprod;                            
00483   int N2[ths->d];
00484 
00485   ths->aN1 = (int*) nfft_malloc(ths->d*sizeof(int));
00486 
00487   ths->a = (double*) nfft_malloc(ths->d*sizeof(double));
00488 
00489   ths->sigma = (double*) nfft_malloc(ths->d*sizeof(double));
00490 
00491   ths->n = ths->N1;
00492 
00493   ths->aN1_total=1;
00494 
00495   for(t = 0; t<ths->d; t++) {
00496     ths->a[t] = 1.0 + (2.0*((double)ths->m))/((double)ths->N1[t]);
00497     ths->aN1[t] = ths->a[t] * ((double)ths->N1[t]);
00498     /* aN1 should be even */
00499     if(ths->aN1[t]%2 != 0)
00500       ths->aN1[t] = ths->aN1[t] +1;
00501 
00502     ths->aN1_total*=ths->aN1[t];
00503     ths->sigma[t] = ((double) ths->N1[t] )/((double) ths->N[t]);;
00504     
00505     /* take the same oversampling factor in the inner NFFT */
00506     N2[t] = ceil(ths->sigma[t]*(ths->aN1[t]));
00507     
00508     /* N2 should be even */
00509     if(N2[t]%2 != 0)
00510       N2[t] = N2[t] +1;
00511   }
00512 
00513   WINDOW_HELP_INIT
00514 
00515   if(ths->nnfft_flags & MALLOC_X)
00516     ths->x = (double*)nfft_malloc(ths->d*ths->M_total*sizeof(double));
00517   if(ths->nnfft_flags & MALLOC_F)
00518     ths->f=(double _Complex*)nfft_malloc(ths->M_total*sizeof(double _Complex));
00519 
00520   if(ths->nnfft_flags & MALLOC_V)
00521     ths->v = (double*)nfft_malloc(ths->d*ths->N_total*sizeof(double));
00522   if(ths->nnfft_flags & MALLOC_F_HAT)
00523     ths->f_hat = (double _Complex*)nfft_malloc(ths->N_total*sizeof(double _Complex));
00524 
00525   if(ths->nnfft_flags & PRE_LIN_PSI)
00526   {
00527     ths->K=(1U<< 10)*(ths->m+1);
00528     ths->psi = (double*) nfft_malloc((ths->K+1)*ths->d*sizeof(double));
00529   }
00530 
00531   if(ths->nnfft_flags & PRE_PSI)
00532     ths->psi = (double*)nfft_malloc(ths->N_total*ths->d*(2*ths->m+2)*sizeof(double));
00533 
00534   if(ths->nnfft_flags & PRE_FULL_PSI)
00535   {
00536       for(t=0,lprod = 1; t<ths->d; t++)
00537           lprod *= 2*ths->m+2;
00538 
00539       ths->psi = (double*)nfft_malloc(ths->N_total*lprod*sizeof(double));
00540 
00541       ths->psi_index_f = (int*) nfft_malloc(ths->N_total*sizeof(int));
00542       ths->psi_index_g = (int*) nfft_malloc(ths->N_total*lprod*sizeof(int));
00543   }
00544 
00545   ths->direct_plan = (nfft_plan*)nfft_malloc(sizeof(nfft_plan));
00546 
00547   nfft_init_guru(ths->direct_plan, ths->d, ths->aN1, ths->M_total, N2, m2,
00548      nfft_flags, fftw_flags);
00549 
00550   ths->direct_plan->x = ths->x;
00551   ths->direct_plan->f = ths->f;
00552   ths->F = ths->direct_plan->f_hat;
00553 
00554   ths->mv_trafo = (void (*) (void* ))nnfft_trafo;
00555   ths->mv_adjoint = (void (*) (void* ))nnfft_adjoint;
00556 }
00557 
00558 void nnfft_init_guru(nnfft_plan *ths, int d, int N_total, int M_total, int *N, int *N1,
00559          int m, unsigned nnfft_flags)
00560 {
00561   int t;                             
00563   unsigned nfft_flags;
00564   unsigned fftw_flags;
00565 
00566   ths->d= d;
00567   ths->M_total= M_total;
00568   ths->N_total= N_total;
00569   ths->m= m;
00570   ths->nnfft_flags= nnfft_flags;
00571   fftw_flags= FFTW_ESTIMATE| FFTW_DESTROY_INPUT;
00572   nfft_flags= PRE_PHI_HUT| MALLOC_F_HAT| FFTW_INIT| FFT_OUT_OF_PLACE;
00573 
00574   if(ths->nnfft_flags & PRE_PSI)
00575     nfft_flags = nfft_flags | PRE_PSI;
00576 
00577   if(ths->nnfft_flags & PRE_FULL_PSI)
00578     nfft_flags = nfft_flags | PRE_FULL_PSI;
00579 
00580   if(ths->nnfft_flags & PRE_LIN_PSI)
00581     nfft_flags = nfft_flags | PRE_LIN_PSI;
00582 
00583   ths->N = (int*) nfft_malloc(ths->d*sizeof(int));
00584   ths->N1 = (int*) nfft_malloc(ths->d*sizeof(int));
00585 
00586   for(t=0; t<d; t++) {
00587     ths->N[t] = N[t];
00588     ths->N1[t] = N1[t];
00589   }
00590   nnfft_init_help(ths,m,nfft_flags,fftw_flags);
00591 }
00592 
00593 void nnfft_init(nnfft_plan *ths, int d, int N_total, int M_total, int *N)
00594 {
00595   int t;                            
00597   unsigned nfft_flags;
00598   unsigned fftw_flags;
00599 
00600   ths->d = d;
00601   ths->M_total = M_total;
00602   ths->N_total = N_total;
00603 
00604   /* m should be greater to get the same accuracy as the nfft */
00605   WINDOW_HELP_ESTIMATE_m;
00606 
00607   ths->N = (int*) nfft_malloc(ths->d*sizeof(int));
00608   ths->N1 = (int*) nfft_malloc(ths->d*sizeof(int));
00609 
00610   for(t=0; t<d; t++) {
00611     ths->N[t] = N[t];
00612 
00613     /* the standard oversampling factor in the nnfft is 1.5 */
00614     ths->N1[t] = ceil(1.5*ths->N[t]);
00615     
00616     /* N1 should be even */
00617     if(ths->N1[t]%2 != 0)
00618       ths->N1[t] = ths->N1[t] +1;
00619   }
00620 
00621   ths->nnfft_flags=PRE_PSI| PRE_PHI_HUT| MALLOC_X| MALLOC_V| MALLOC_F_HAT| MALLOC_F;
00622   nfft_flags= PRE_PSI| PRE_PHI_HUT| MALLOC_F_HAT| FFTW_INIT| FFT_OUT_OF_PLACE;
00623 
00624   fftw_flags= FFTW_ESTIMATE| FFTW_DESTROY_INPUT;
00625 
00626   nnfft_init_help(ths,ths->m,nfft_flags,fftw_flags);
00627 }
00628 
00629 void nnfft_finalize(nnfft_plan *ths)
00630 {
00631 
00632   nfft_finalize(ths->direct_plan);
00633 
00634   nfft_free(ths->direct_plan);
00635 
00636   nfft_free(ths->aN1);
00637   nfft_free(ths->N);
00638   nfft_free(ths->N1);
00639 
00640   if(ths->nnfft_flags & PRE_FULL_PSI)
00641     {
00642       nfft_free(ths->psi_index_g);
00643       nfft_free(ths->psi_index_f);
00644       nfft_free(ths->psi);
00645     }
00646 
00647   if(ths->nnfft_flags & PRE_PSI)
00648     nfft_free(ths->psi);
00649 
00650   if(ths->nnfft_flags & PRE_LIN_PSI)
00651     nfft_free(ths->psi);
00652 
00653   if(ths->nnfft_flags & PRE_PHI_HUT)
00654     nfft_free(ths->c_phi_inv);
00655 
00656   if(ths->nnfft_flags & MALLOC_F)
00657     nfft_free(ths->f);
00658 
00659   if(ths->nnfft_flags & MALLOC_F_HAT)
00660     nfft_free(ths->f_hat);
00661 
00662   if(ths->nnfft_flags & MALLOC_X)
00663     nfft_free(ths->x);
00664 
00665   if(ths->nnfft_flags & MALLOC_V)
00666     nfft_free(ths->v);
00667 }

Generated on 17 Aug 2009 by Doxygen 1.5.3