NFFT Logo 3.1.0 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 3100 2009-03-12 08:42:48Z 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 "util.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*               \
00035                                               sizeof(double _Complex));
00036 #define MACRO_nndft_init_result_transposed MACRO_nndft_init_result_adjoint
00037 
00038 #define MACRO_nndft_sign_trafo      (-2.0*PI)
00039 #define MACRO_nndft_sign_conjugated (+2.0*PI)
00040 #define MACRO_nndft_sign_adjoint    (+2.0*PI)
00041 #define MACRO_nndft_sign_transposed (-2.0*PI)
00042 
00043 #define MACRO_nndft_compute_trafo (*fj) += (*f_hat_k)*cexp(+I*omega);
00044 
00045 #define MACRO_nndft_compute_conjugated MACRO_nndft_compute_trafo
00046 
00047 #define MACRO_nndft_compute_adjoint (*f_hat_k) += (*fj)*cexp(+ _Complex_I*omega);
00048 
00049 #define MACRO_nndft_compute_transposed MACRO_nndft_compute_adjoint
00050 
00051 #define MACRO_nndft(which_one)                                                 \
00052 void nndft_ ## which_one (nnfft_plan *ths)                                    \
00053 {                                                                              \
00054   int j;                                \
00055   int t;                                \
00056   int l;                                \
00057   double _Complex *f_hat, *f;                   \
00058   double _Complex *f_hat_k;                     \
00059   double _Complex *fj;                          \
00060   double omega;                         \
00061                                                                                \
00062   f_hat=ths->f_hat; f=ths->f;                                                \
00063                                                                                \
00064   MACRO_nndft_init_result_ ## which_one                                        \
00065                                                                                \
00066   for(j=0, fj=f; j<ths->M_total; j++, fj++)                                     \
00067   {                                                                            \
00068     for(l=0, f_hat_k=f_hat; l<ths->N_total; l++, f_hat_k++)                    \
00069     {                                                                          \
00070       omega=0.0;                                                               \
00071       for(t = 0; t<ths->d; t++)                                               \
00072         omega+=ths->v[l*ths->d+t] * ths->x[j*ths->d+t] * ths->N[t];           \
00073                                                                                \
00074       omega*= MACRO_nndft_sign_ ## which_one;                                  \
00075                                                                                \
00076       MACRO_nndft_compute_ ## which_one                                        \
00077                                                                                \
00078      } /* for(l) */                                                            \
00079    } /* for(j) */                                                              \
00080 } /* nndft_trafo */                                                            \
00081 
00082 MACRO_nndft(trafo)
00083 MACRO_nndft(adjoint)
00084 
00087 void nnfft_uo(nnfft_plan *ths,int j,int *up,int *op,int act_dim)
00088 {
00089   double c;
00090   int u,o;
00091 
00092   c = ths->v[j*ths->d+act_dim] * ths->n[act_dim];
00093 
00094   u = c; o = c;
00095   if(c < 0)
00096     u = u-1;
00097   else
00098     o = o+1;
00099 
00100   u = u - (ths->m); o = o + (ths->m);
00101 
00102   up[0]=u; op[0]=o;
00103 }
00104 
00108 #define MACRO_nnfft_B_init_result_A memset(f,0,ths->M_total*sizeof(double _Complex));
00109 #define MACRO_nnfft_B_init_result_T memset(g,0,ths->aN1_total*sizeof(double _Complex));
00110 
00111 #define MACRO_nnfft_B_PRE_FULL_PSI_compute_A {                                  \
00112   (*fj) += ths->psi[ix] * g[ths->psi_index_g[ix]];                           \
00113 }
00114 
00115 #define MACRO_nnfft_B_PRE_FULL_PSI_compute_T {                                  \
00116   g[ths->psi_index_g[ix]] += ths->psi[ix] * (*fj);                           \
00117 }
00118 
00119 #define MACRO_nnfft_B_compute_A {                                               \
00120   (*fj) += phi_prod[ths->d] * g[ll_plain[ths->d]];                        \
00121 }
00122 
00123 #define MACRO_nnfft_B_compute_T {                                               \
00124   g[ll_plain[ths->d]] += phi_prod[ths->d] * (*fj);                        \
00125 }
00126 
00127   /* Gewicht, d.h., Nachkommaanteil y-y_u im Speicher halten!!! */
00128 #define MACRO_with_PRE_LIN_PSI (ths->psi[(ths->K+1)*t2+y_u[t2]]*             \
00129                                 (y_u[t2]+1-y[t2]) +                            \
00130                                 ths->psi[(ths->K+1)*t2+y_u[t2]+1]*           \
00131                                 (y[t2]-y_u[t2]))
00132 #define MACRO_with_PRE_PSI     ths->psi[(j*ths->d+t2)*(2*ths->m+2)+lj[t2]]
00133 #define MACRO_without_PRE_PSI  PHI(-ths->v[j*ths->d+t2]+                      \
00134                                ((double)l[t2])/ths->N1[t2], t2)
00135 
00136 #define MACRO_init_uo_l_lj_t {                                                 \
00137   for(t = ths->d-1; t>=0; t--)                                                \
00138     {                                                                          \
00139       nnfft_uo(ths,j,&u[t],&o[t],t);                                           \
00140       l[t] = u[t];                                                             \
00141       lj[t] = 0;                                                               \
00142     } /* for(t) */                                                             \
00143   t++;                                                                         \
00144 }
00145 
00146 #define MACRO_update_with_PRE_PSI_LIN {                                        \
00147   for(t2=t; t2<ths->d; t2++)                                                  \
00148     {                                                                          \
00149       y[t2] = fabs(((-ths->N1[t2]*ths->v[j*ths->d+t2]+(double)l[t2])          \
00150           * ((double)ths->K))/(ths->m+1));                                   \
00151       y_u[t2] = (int)y[t2];                                                    \
00152     } /* for(t2) */                                                            \
00153 }
00154 
00155 #define MACRO_update_phi_prod_ll_plain(which_one) {                            \
00156   for(t2=t; t2<ths->d; t2++)                                                  \
00157     {                                                                          \
00158       phi_prod[t2+1]=phi_prod[t2]* MACRO_ ## which_one;                        \
00159       ll_plain[t2+1]=ll_plain[t2]*ths->aN1[t2] +(l[t2]+ths->aN1[t2]*3/2)%ths->aN1[t2]; /* 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     for(j=0; j<ths->M_total; j++)
00273       ths->f[j] *= ths->c_phi_inv[j];
00274   } else {
00275     for(j=0; j<ths->M_total; j++)
00276     {
00277       tmp = 1.0;
00278       /* multiply with N1, because x was modified */
00279       for(t=0; t<ths->d; t++)
00280         tmp*= 1.0 /((PHI_HUT(ths->x[ths->d*j + t]*((double)ths->N[t]),t)) );
00281       ths->f[j] *= tmp;
00282     }
00283   }
00284 }
00285 
00288 void nnfft_trafo(nnfft_plan *ths)
00289 {
00290   int j,t;
00291 
00292   nnfft_B_T(ths);
00293 
00294   for(j=0;j<ths->M_total;j++) {
00295     for(t=0;t<ths->d;t++) {
00296       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
00297     }
00298   }
00299 
00300   ths->direct_plan->f = ths->f;
00301   nfft_trafo(ths->direct_plan);
00302 
00303   for(j=0;j<ths->M_total;j++) {
00304     for(t=0;t<ths->d;t++) {
00305       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
00306     }
00307   }
00308 
00309   nnfft_D(ths);
00310 } /* nnfft_trafo */
00311 
00312 void nnfft_adjoint(nnfft_plan *ths)
00313 {
00314   int j,t;
00315 
00316   nnfft_D(ths);
00317 
00318   for(j=0;j<ths->M_total;j++) {
00319     for(t=0;t<ths->d;t++) {
00320       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
00321     }
00322   }
00323 
00324   ths->direct_plan->f = ths->f;
00325   nfft_adjoint(ths->direct_plan);
00326 
00327   for(j=0;j<ths->M_total;j++) {
00328     for(t=0;t<ths->d;t++) {
00329       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
00330     }
00331   }
00332 
00333   nnfft_B_A(ths);
00334 } /* nnfft_adjoint */
00335 
00338 void nnfft_precompute_phi_hut(nnfft_plan *ths)
00339 {
00340   int j;                                
00341   int t;                                
00342   double tmp;
00343 
00344   ths->c_phi_inv= (double*)nfft_malloc(ths->M_total*sizeof(double));
00345 
00346   for(j=0; j<ths->M_total; j++)
00347     {
00348       tmp = 1.0;
00349       for(t=0; t<ths->d; t++)
00350         tmp*= 1.0 /(PHI_HUT(ths->x[ths->d*j + t]*((double)ths->N[t]),t));
00351       ths->c_phi_inv[j]=tmp;
00352     }
00353 } /* nnfft_phi_hut */
00354 
00355 
00360 void nnfft_precompute_lin_psi(nnfft_plan *ths)
00361 {
00362   int t;                                
00363   int j;                                
00364   double step;                          
00366   nfft_precompute_lin_psi(ths->direct_plan);
00367 
00368   for (t=0; t<ths->d; t++)
00369     {
00370       step=((double)(ths->m+1))/(ths->K*ths->N1[t]);
00371       for(j=0;j<=ths->K;j++)
00372         {
00373           ths->psi[(ths->K+1)*t + j] = PHI(j*step,t);
00374         } /* for(j) */
00375     } /* for(t) */
00376 }
00377 
00378 void nnfft_precompute_psi(nnfft_plan *ths)
00379 {
00380   int t;                                
00381   int j;                                
00382   int l;                                
00383   int lj;                               
00384   int u, o;                             
00386   for (t=0; t<ths->d; t++)
00387     for(j=0;j<ths->N_total;j++)
00388       {
00389         nnfft_uo(ths,j,&u,&o,t);
00390 
00391         for(l=u, lj=0; l <= o; l++, lj++)
00392           ths->psi[(j*ths->d+t)*(2*ths->m+2)+lj]=
00393             (PHI((-ths->v[j*ths->d+t]+((double)l)/((double)ths->N1[t])),t));
00394       } /* for(j) */
00395 
00396   for(j=0;j<ths->M_total;j++) {
00397     for(t=0;t<ths->d;t++) {
00398       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
00399     }
00400   }
00401 
00402   nfft_precompute_psi(ths->direct_plan);
00403 
00404   for(j=0;j<ths->M_total;j++) {
00405     for(t=0;t<ths->d;t++) {
00406       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
00407     }
00408   }
00409   /* for(t) */
00410 } /* nfft_precompute_psi */
00411 
00412 
00413 
00417 void nnfft_precompute_full_psi(nnfft_plan *ths)
00418 {
00419   int t,t2;                             
00420   int j;                                
00421   int l_L;                              
00422   int l[ths->d];                       
00423   int lj[ths->d];                      
00424   int ll_plain[ths->d+1];              
00425   int lprod;                            
00426   int u[ths->d], o[ths->d];           
00428   double phi_prod[ths->d+1];
00429 
00430   int ix,ix_old;
00431 
00432   for(j=0;j<ths->M_total;j++) {
00433     for(t=0;t<ths->d;t++) {
00434       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
00435     }
00436   }
00437 
00438   nnfft_precompute_psi(ths);
00439 
00440   nfft_precompute_full_psi(ths->direct_plan);
00441 
00442   for(j=0;j<ths->M_total;j++) {
00443     for(t=0;t<ths->d;t++) {
00444       ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
00445     }
00446   }
00447 
00448   phi_prod[0]=1;
00449   ll_plain[0]=0;
00450 
00451   for(t=0,lprod = 1; t<ths->d; t++)
00452     lprod *= 2*ths->m+2;
00453 
00454   for(j=0,ix=0,ix_old=0; j<ths->N_total; j++)
00455     {
00456       MACRO_init_uo_l_lj_t;
00457 
00458       for(l_L=0; l_L<lprod; l_L++, ix++)
00459         {
00460           MACRO_update_phi_prod_ll_plain(without_PRE_PSI);
00461 
00462           ths->psi_index_g[ix]=ll_plain[ths->d];
00463           ths->psi[ix]=phi_prod[ths->d];
00464 
00465           MACRO_count_uo_l_lj_t;
00466         } /* for(l_L) */
00467 
00468 
00469       ths->psi_index_f[j]=ix-ix_old;
00470       ix_old=ix;
00471     } /* for(j) */
00472 }
00473 
00474 void nnfft_init_help(nnfft_plan *ths, int m2, unsigned nfft_flags, unsigned fftw_flags)
00475 {
00476   int t;                                
00477   int lprod;                            
00478   int N2[ths->d];
00479 
00480   ths->aN1 = (int*) nfft_malloc(ths->d*sizeof(int));
00481 
00482   ths->a = (double*) nfft_malloc(ths->d*sizeof(double));
00483 
00484   ths->sigma = (double*) nfft_malloc(ths->d*sizeof(double));
00485 
00486   ths->n = ths->N1;
00487 
00488   ths->aN1_total=1;
00489 
00490   for(t = 0; t<ths->d; t++) {
00491     ths->a[t] = 1.0 + (2.0*((double)ths->m))/((double)ths->N1[t]);
00492     ths->aN1[t] = ths->a[t] * ((double)ths->N1[t]);
00493     /* aN1 should be even */
00494     if(ths->aN1[t]%2 != 0)
00495       ths->aN1[t] = ths->aN1[t] +1;
00496 
00497     ths->aN1_total*=ths->aN1[t];
00498     ths->sigma[t] = ((double) ths->N1[t] )/((double) ths->N[t]);;
00499 
00500     /* take the same oversampling factor in the inner NFFT */
00501     N2[t] = ceil(ths->sigma[t]*(ths->aN1[t]));
00502 
00503     /* N2 should be even */
00504     if(N2[t]%2 != 0)
00505       N2[t] = N2[t] +1;
00506   }
00507 
00508   WINDOW_HELP_INIT
00509 
00510   if(ths->nnfft_flags & MALLOC_X)
00511     ths->x = (double*)nfft_malloc(ths->d*ths->M_total*
00512                                         sizeof(double));
00513 
00514   if(ths->nnfft_flags & MALLOC_V)
00515     ths->v = (double*)nfft_malloc(ths->d*ths->N_total*
00516                                         sizeof(double));
00517 
00518   if(ths->nnfft_flags & MALLOC_F_HAT)
00519     ths->f_hat = (double _Complex*)nfft_malloc(ths->N_total*
00520                                                   sizeof(double _Complex));
00521   if(ths->nnfft_flags & MALLOC_F)
00522     ths->f=(double _Complex*)nfft_malloc(ths->M_total*sizeof(double _Complex));
00523 
00524   if(ths->nnfft_flags & PRE_LIN_PSI)
00525   {
00526     ths->K=100000; /* estimate is badly needed !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!*/
00527     ths->psi = (double*) nfft_malloc((ths->K+1)*ths->d*sizeof(double));
00528   }
00529 
00530   /* NO FFTW_MALLOC HERE */
00531   if(ths->nnfft_flags & PRE_PSI)
00532     ths->psi = (double*)nfft_malloc(ths->N_total*ths->d*
00533                                            (2*ths->m+2)*sizeof(double));
00534 
00535   if(ths->nnfft_flags & PRE_FULL_PSI)
00536   {
00537       for(t=0,lprod = 1; t<ths->d; t++)
00538           lprod *= 2*ths->m+2;
00539 
00540       ths->psi = (double*) nfft_malloc(ths->M_total*lprod*sizeof(double));
00541 
00542       ths->psi_index_f = (int*) nfft_malloc(ths->M_total*sizeof(int));
00543       ths->psi_index_g = (int*) nfft_malloc(ths->M_total*lprod*sizeof(int));
00544   }
00545 
00546   ths->direct_plan = (nfft_plan*)nfft_malloc(sizeof(nfft_plan));
00547 
00548   nfft_init_guru(ths->direct_plan, ths->d, ths->aN1, ths->M_total, N2, m2,
00549                  nfft_flags, fftw_flags);
00550 
00551   ths->direct_plan->x = ths->x;
00552   ths->direct_plan->f = ths->f;
00553   ths->F = ths->direct_plan->f_hat;
00554 
00555   ths->mv_trafo = (void (*) (void* ))nnfft_trafo;
00556   ths->mv_adjoint = (void (*) (void* ))nnfft_adjoint;
00557 }
00558 
00559 void nnfft_init_guru(nnfft_plan *ths, int d, int N_total, int M_total, int *N, int *N1,
00560                         int m, unsigned nnfft_flags)
00561 {
00562   int t;                             
00564   unsigned nfft_flags;
00565   unsigned fftw_flags;
00566 
00567   ths->d= d;
00568   ths->M_total= M_total;
00569   ths->N_total= N_total;
00570   ths->m= m;
00571   ths->nnfft_flags= nnfft_flags;
00572   fftw_flags= FFTW_ESTIMATE| FFTW_DESTROY_INPUT;
00573   nfft_flags= PRE_PHI_HUT| MALLOC_F_HAT| FFTW_INIT| FFT_OUT_OF_PLACE;
00574 
00575   if(ths->nnfft_flags & PRE_PSI)
00576     nfft_flags = nfft_flags | PRE_PSI;
00577 
00578   if(ths->nnfft_flags & PRE_FULL_PSI)
00579     nfft_flags = nfft_flags | PRE_FULL_PSI;
00580 
00581   if(ths->nnfft_flags & PRE_LIN_PSI)
00582     nfft_flags = nfft_flags | PRE_LIN_PSI;
00583 
00584   ths->N = (int*) nfft_malloc(ths->d*sizeof(int));
00585   ths->N1 = (int*) nfft_malloc(ths->d*sizeof(int));
00586 
00587   for(t=0; t<d; t++) {
00588     ths->N[t] = N[t];
00589     ths->N1[t] = N1[t];
00590   }
00591   nnfft_init_help(ths,m,nfft_flags,fftw_flags);
00592 }
00593 
00594 void nnfft_init(nnfft_plan *ths, int d, int N_total, int M_total, int *N)
00595 {
00596   int t;                            
00598   unsigned nfft_flags;
00599   unsigned fftw_flags;
00600 
00601   ths->d = d;
00602   ths->M_total = M_total;
00603   ths->N_total = N_total;
00604 
00605   /* m should be greater to get the same accuracy as the nfft */
00606   WINDOW_HELP_ESTIMATE_m;
00607 
00608   ths->N = (int*) nfft_malloc(ths->d*sizeof(int));
00609   ths->N1 = (int*) nfft_malloc(ths->d*sizeof(int));
00610 
00611   for(t=0; t<d; t++) {
00612     ths->N[t] = N[t];
00613 
00614     /* the standard oversampling factor in the nnfft is 1.5 */
00615     ths->N1[t] = ceil(1.5*ths->N[t]);
00616 
00617     /* N1 should be even */
00618     if(ths->N1[t]%2 != 0)
00619       ths->N1[t] = ths->N1[t] +1;
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   nfft_finalize(ths->direct_plan);
00632 
00633   nfft_free(ths->direct_plan);
00634 
00635   nfft_free(ths->aN1);
00636   nfft_free(ths->N);
00637   nfft_free(ths->N1);
00638 
00639   if(ths->nnfft_flags & PRE_FULL_PSI)
00640     {
00641       nfft_free(ths->psi_index_g);
00642       nfft_free(ths->psi_index_f);
00643       nfft_free(ths->psi);
00644     }
00645 
00646   if(ths->nnfft_flags & PRE_PSI)
00647     nfft_free(ths->psi);
00648 
00649   if(ths->nnfft_flags & PRE_LIN_PSI)
00650     nfft_free(ths->psi);
00651 
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 19 Mar 2009 by Doxygen 1.5.3