NFFT Logo 3.2.3
nnfft.c
1 /*
2  * Copyright (c) 2002, 2012 Jens Keiner, Stefan Kunis, Daniel Potts
3  *
4  * This program is free software; you can redistribute it and/or modify it under
5  * the terms of the GNU General Public License as published by the Free Software
6  * Foundation; either version 2 of the License, or (at your option) any later
7  * version.
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
12  * details.
13  *
14  * You should have received a copy of the GNU General Public License along with
15  * this program; if not, write to the Free Software Foundation, Inc., 51
16  * Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
17  */
18 
19 /* $Id: nnfft.c 3896 2012-10-10 12:19:26Z tovo $ */
20 
21 #include "config.h"
22 
23 #include <stdio.h>
24 #include <math.h>
25 #include <string.h>
26 #include <stdlib.h>
27 #ifdef HAVE_COMPLEX_H
28 #include <complex.h>
29 #endif
30 #include "nfft3util.h"
31 #include "nfft3.h"
32 #include "infft.h"
33 
34 
35 #define MACRO_nndft_init_result_trafo memset(f,0,ths->M_total*sizeof(double _Complex));
36 #define MACRO_nndft_init_result_conjugated MACRO_nndft_init_result_trafo
37 #define MACRO_nndft_init_result_adjoint memset(f_hat,0,ths->N_total*sizeof(double _Complex));
38 #define MACRO_nndft_init_result_transposed MACRO_nndft_init_result_adjoint
39 
40 #define MACRO_nndft_sign_trafo (-2.0*PI)
41 #define MACRO_nndft_sign_conjugated (+2.0*PI)
42 #define MACRO_nndft_sign_adjoint (+2.0*PI)
43 #define MACRO_nndft_sign_transposed (-2.0*PI)
44 
45 #define MACRO_nndft_compute_trafo (*fj) += (*f_hat_k)*cexp(+ _Complex_I*omega);
46 
47 #define MACRO_nndft_compute_conjugated MACRO_nndft_compute_trafo
48 
49 #define MACRO_nndft_compute_adjoint (*f_hat_k) += (*fj)*cexp(+ _Complex_I*omega);
50 
51 #define MACRO_nndft_compute_transposed MACRO_nndft_compute_adjoint
52 
53 #define MACRO_nndft(which_one) \
54 void nnfft_ ## which_one ## _direct (nnfft_plan *ths) \
55 { \
56  int j; \
57  int t; \
58  int l; \
59  double _Complex *f_hat, *f; \
60  double _Complex *f_hat_k; \
61  double _Complex *fj; \
62  double omega; \
63  \
64  f_hat=ths->f_hat; f=ths->f; \
65  \
66  MACRO_nndft_init_result_ ## which_one \
67  \
68  for(j=0, fj=f; j<ths->M_total; j++, fj++) \
69  { \
70  for(l=0, f_hat_k=f_hat; l<ths->N_total; l++, f_hat_k++) \
71  { \
72  omega=0.0; \
73  for(t = 0; t<ths->d; t++) \
74  omega+=ths->v[l*ths->d+t] * ths->x[j*ths->d+t] * ths->N[t]; \
75  \
76  omega*= MACRO_nndft_sign_ ## which_one; \
77  \
78  MACRO_nndft_compute_ ## which_one \
79  \
80  } /* for(l) */ \
81  } /* for(j) */ \
82 } /* nndft_trafo */ \
83 
84 MACRO_nndft(trafo)
85 MACRO_nndft(adjoint)
86 
89 static void nnfft_uo(nnfft_plan *ths,int j,int *up,int *op,int act_dim)
90 {
91  double c;
92  int u,o;
93 
94  c = ths->v[j*ths->d+act_dim] * ths->n[act_dim];
95 
96  u = c; o = c;
97  if(c < 0)
98  u = u-1;
99  else
100  o = o+1;
101 
102  u = u - (ths->m); o = o + (ths->m);
103 
104  up[0]=u; op[0]=o;
105 }
106 
110 #define MACRO_nnfft_B_init_result_A memset(f,0,ths->N_total*sizeof(double _Complex));
111 #define MACRO_nnfft_B_init_result_T memset(g,0,ths->aN1_total*sizeof(double _Complex));
112 
113 #define MACRO_nnfft_B_PRE_FULL_PSI_compute_A { \
114  (*fj) += ths->psi[ix] * g[ths->psi_index_g[ix]]; \
115 }
116 
117 #define MACRO_nnfft_B_PRE_FULL_PSI_compute_T { \
118  g[ths->psi_index_g[ix]] += ths->psi[ix] * (*fj); \
119 }
120 
121 #define MACRO_nnfft_B_compute_A { \
122  (*fj) += phi_prod[ths->d] * g[ll_plain[ths->d]]; \
123 }
124 
125 #define MACRO_nnfft_B_compute_T { \
126  g[ll_plain[ths->d]] += phi_prod[ths->d] * (*fj); \
127 }
128 
129 #define MACRO_with_PRE_LIN_PSI (ths->psi[(ths->K+1)*t2+y_u[t2]]* \
130  (y_u[t2]+1-y[t2]) + \
131  ths->psi[(ths->K+1)*t2+y_u[t2]+1]* \
132  (y[t2]-y_u[t2]))
133 #define MACRO_with_PRE_PSI ths->psi[(j*ths->d+t2)*(2*ths->m+2)+lj[t2]]
134 #define MACRO_without_PRE_PSI PHI(-ths->v[j*ths->d+t2]+ \
135  ((double)l[t2])/ths->N1[t2], t2)
136 
137 #define MACRO_init_uo_l_lj_t { \
138  for(t = ths->d-1; t>=0; t--) \
139  { \
140  nnfft_uo(ths,j,&u[t],&o[t],t); \
141  l[t] = u[t]; \
142  lj[t] = 0; \
143  } /* for(t) */ \
144  t++; \
145 }
146 
147 #define MACRO_update_with_PRE_PSI_LIN { \
148  for(t2=t; t2<ths->d; t2++) \
149  { \
150  y[t2] = fabs(((-ths->N1[t2]*ths->v[j*ths->d+t2]+(double)l[t2]) \
151  * ((double)ths->K))/(ths->m+1)); \
152  y_u[t2] = (int)y[t2]; \
153  } /* for(t2) */ \
154 }
155 
156 #define MACRO_update_phi_prod_ll_plain(which_one) { \
157  for(t2=t; t2<ths->d; t2++) \
158  { \
159  phi_prod[t2+1]=phi_prod[t2]* MACRO_ ## which_one; \
160  ll_plain[t2+1]=ll_plain[t2]*ths->aN1[t2] + \
161  (l[t2]+ths->aN1[t2]*3/2)%ths->aN1[t2]; \
162  /* 3/2 because of the (not needed) fftshift and to be in [0 aN1[t2]]?!*/\
163  } /* for(t2) */ \
164 }
165 
166 #define MACRO_count_uo_l_lj_t { \
167  for(t = ths->d-1; (t>0)&&(l[t]==o[t]); t--) \
168  { \
169  l[t] = u[t]; \
170  lj[t] = 0; \
171  } /* for(t) */ \
172  \
173  l[t]++; \
174  lj[t]++; \
175 }
176 
177 #define MACRO_nnfft_B(which_one) \
178 static inline void nnfft_B_ ## which_one (nnfft_plan *ths) \
179 { \
180  int lprod; \
181  int u[ths->d], o[ths->d]; \
182  int t, t2; \
183  int j; \
184  int l_L, ix; \
185  int l[ths->d]; \
186  int lj[ths->d]; \
187  int ll_plain[ths->d+1]; \
188  double phi_prod[ths->d+1]; \
189  double _Complex *f, *g; \
190  double _Complex *fj; \
191  double y[ths->d]; \
192  int y_u[ths->d]; \
193  \
194  f=ths->f_hat; g=ths->F; \
195  \
196  MACRO_nnfft_B_init_result_ ## which_one \
197  \
198  if(ths->nnfft_flags & PRE_FULL_PSI) \
199  { \
200  for(ix=0, j=0, fj=f; j<ths->N_total; j++,fj++) \
201  for(l_L=0; l_L<ths->psi_index_f[j]; l_L++, ix++) \
202  MACRO_nnfft_B_PRE_FULL_PSI_compute_ ## which_one; \
203  return; \
204  } \
205  \
206  phi_prod[0]=1; \
207  ll_plain[0]=0; \
208  \
209  for(t=0,lprod = 1; t<ths->d; t++) \
210  lprod *= (2*ths->m+2); \
211  \
212  if(ths->nnfft_flags & PRE_PSI) \
213  { \
214  for(j=0, fj=f; j<ths->N_total; j++, fj++) \
215  { \
216  MACRO_init_uo_l_lj_t; \
217  \
218  for(l_L=0; l_L<lprod; l_L++) \
219  { \
220  MACRO_update_phi_prod_ll_plain(with_PRE_PSI); \
221  \
222  MACRO_nnfft_B_compute_ ## which_one; \
223  \
224  MACRO_count_uo_l_lj_t; \
225  } /* for(l_L) */ \
226  } /* for(j) */ \
227  return; \
228  } /* if(PRE_PSI) */ \
229  \
230  if(ths->nnfft_flags & PRE_LIN_PSI) \
231  { \
232  for(j=0, fj=f; j<ths->N_total; j++, fj++) \
233  { \
234  MACRO_init_uo_l_lj_t; \
235  \
236  for(l_L=0; l_L<lprod; l_L++) \
237  { \
238  MACRO_update_with_PRE_PSI_LIN; \
239  \
240  MACRO_update_phi_prod_ll_plain(with_PRE_LIN_PSI); \
241  \
242  MACRO_nnfft_B_compute_ ## which_one; \
243  \
244  MACRO_count_uo_l_lj_t; \
245  } /* for(l_L) */ \
246  } /* for(j) */ \
247  return; \
248  } /* if(PRE_LIN_PSI) */ \
249  \
250  /* no precomputed psi at all */ \
251  for(j=0, fj=f; j<ths->N_total; j++, fj++) \
252  { \
253  \
254  MACRO_init_uo_l_lj_t; \
255  \
256  for(l_L=0; l_L<lprod; l_L++) \
257  { \
258  MACRO_update_phi_prod_ll_plain(without_PRE_PSI); \
259  \
260  MACRO_nnfft_B_compute_ ## which_one; \
261  \
262  MACRO_count_uo_l_lj_t; \
263  } /* for(l_L) */ \
264  } /* for(j) */ \
265 } /* nnfft_B */
266 
267 MACRO_nnfft_B(A)
268 MACRO_nnfft_B(T)
269 
270 static inline void nnfft_D (nnfft_plan *ths){
271  int j,t;
272  double tmp;
273 
274  if(ths->nnfft_flags & PRE_PHI_HUT)
275  {
276  for(j=0; j<ths->M_total; j++)
277  ths->f[j] *= ths->c_phi_inv[j];
278  }
279  else
280  {
281  for(j=0; j<ths->M_total; j++)
282  {
283  tmp = 1.0;
284  /* multiply with N1, because x was modified */
285  for(t=0; t<ths->d; t++)
286  tmp*= 1.0 /((PHI_HUT(ths->x[ths->d*j + t]*((double)ths->N[t]),t)) );
287  ths->f[j] *= tmp;
288  }
289  }
290 }
291 
295 {
296  int j,t;
297 
298  nnfft_B_T(ths);
299 
300  for(j=0;j<ths->M_total;j++) {
301  for(t=0;t<ths->d;t++) {
302  ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
303  }
304  }
305 
306  /* allows for external swaps of ths->f */
307  ths->direct_plan->f = ths->f;
308 
309  nfft_trafo(ths->direct_plan);
310 
311  for(j=0;j<ths->M_total;j++) {
312  for(t=0;t<ths->d;t++) {
313  ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
314  }
315  }
316 
317  nnfft_D(ths);
318 } /* nnfft_trafo */
319 
320 void nnfft_adjoint(nnfft_plan *ths)
321 {
322  int j,t;
323 
324  nnfft_D(ths);
325 
326  for(j=0;j<ths->M_total;j++) {
327  for(t=0;t<ths->d;t++) {
328  ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
329  }
330  }
331 
332  /* allows for external swaps of ths->f */
333  ths->direct_plan->f=ths->f;
334 
335  nfft_adjoint(ths->direct_plan);
336 
337  for(j=0;j<ths->M_total;j++) {
338  for(t=0;t<ths->d;t++) {
339  ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
340  }
341  }
342 
343  nnfft_B_A(ths);
344 } /* nnfft_adjoint */
345 
349 {
350  int j;
351  int t;
352  double tmp;
353 
354  ths->c_phi_inv= (double*)nfft_malloc(ths->M_total*sizeof(double));
355 
356  for(j=0; j<ths->M_total; j++)
357  {
358  tmp = 1.0;
359  for(t=0; t<ths->d; t++)
360  tmp*= 1.0 /(PHI_HUT(ths->x[ths->d*j + t]*((double)ths->N[t]),t));
361  ths->c_phi_inv[j]=tmp;
362  }
363 } /* nnfft_phi_hut */
364 
365 
369 {
370  int t;
371  int j;
372  double step;
375 
376  for (t=0; t<ths->d; t++)
377  {
378  step=((double)(ths->m+1))/(ths->K*ths->N1[t]);
379  for(j=0;j<=ths->K;j++)
380  {
381  ths->psi[(ths->K+1)*t + j] = PHI(j*step,t);
382  } /* for(j) */
383  } /* for(t) */
384 }
385 
386 void nnfft_precompute_psi(nnfft_plan *ths)
387 {
388  int t;
389  int j;
390  int l;
391  int lj;
392  int u, o;
394  for (t=0; t<ths->d; t++)
395  for(j=0;j<ths->N_total;j++)
396  {
397  nnfft_uo(ths,j,&u,&o,t);
398 
399  for(l=u, lj=0; l <= o; l++, lj++)
400  ths->psi[(j*ths->d+t)*(2*ths->m+2)+lj]=
401  (PHI((-ths->v[j*ths->d+t]+((double)l)/((double)ths->N1[t])),t));
402  } /* for(j) */
403 
404  for(j=0;j<ths->M_total;j++) {
405  for(t=0;t<ths->d;t++) {
406  ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
407  }
408  }
409 
410  nfft_precompute_psi(ths->direct_plan);
411 
412  for(j=0;j<ths->M_total;j++) {
413  for(t=0;t<ths->d;t++) {
414  ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
415  }
416  }
417  /* for(t) */
418 } /* nfft_precompute_psi */
419 
420 
421 
426 {
427  int t,t2;
428  int j;
429  int l_L;
430  int l[ths->d];
431  int lj[ths->d];
432  int ll_plain[ths->d+1];
433  int lprod;
434  int u[ths->d], o[ths->d];
436  double phi_prod[ths->d+1];
437 
438  int ix,ix_old;
439 
440  for(j=0;j<ths->M_total;j++) {
441  for(t=0;t<ths->d;t++) {
442  ths->x[j*ths->d+t]= ths->x[j*ths->d+t] / ((double)ths->sigma[t]);
443  }
444  }
445 
446  nnfft_precompute_psi(ths);
447 
448  nfft_precompute_full_psi(ths->direct_plan);
449 
450  for(j=0;j<ths->M_total;j++) {
451  for(t=0;t<ths->d;t++) {
452  ths->x[j*ths->d+t]= ths->x[j*ths->d+t] * ((double)ths->sigma[t]);
453  }
454  }
455 
456  phi_prod[0]=1;
457  ll_plain[0]=0;
458 
459  for(t=0,lprod = 1; t<ths->d; t++)
460  lprod *= 2*ths->m+2;
461 
462  for(j=0,ix=0,ix_old=0; j<ths->N_total; j++)
463  {
464  MACRO_init_uo_l_lj_t;
465 
466  for(l_L=0; l_L<lprod; l_L++, ix++)
467  {
468  MACRO_update_phi_prod_ll_plain(without_PRE_PSI);
469 
470  ths->psi_index_g[ix]=ll_plain[ths->d];
471  ths->psi[ix]=phi_prod[ths->d];
472 
473  MACRO_count_uo_l_lj_t;
474  } /* for(l_L) */
475 
476 
477  ths->psi_index_f[j]=ix-ix_old;
478  ix_old=ix;
479  } /* for(j) */
480 }
481 
482 static void nnfft_init_help(nnfft_plan *ths, int m2, unsigned nfft_flags, unsigned fftw_flags)
483 {
484  int t;
485  int lprod;
486  int N2[ths->d];
487 
488  ths->aN1 = (int*) nfft_malloc(ths->d*sizeof(int));
489 
490  ths->a = (double*) nfft_malloc(ths->d*sizeof(double));
491 
492  ths->sigma = (double*) nfft_malloc(ths->d*sizeof(double));
493 
494  ths->n = ths->N1;
495 
496  ths->aN1_total=1;
497 
498  for(t = 0; t<ths->d; t++) {
499  ths->a[t] = 1.0 + (2.0*((double)ths->m))/((double)ths->N1[t]);
500  ths->aN1[t] = ths->a[t] * ((double)ths->N1[t]);
501  /* aN1 should be even */
502  if(ths->aN1[t]%2 != 0)
503  ths->aN1[t] = ths->aN1[t] +1;
504 
505  ths->aN1_total*=ths->aN1[t];
506  ths->sigma[t] = ((double) ths->N1[t] )/((double) ths->N[t]);;
507 
508  /* take the same oversampling factor in the inner NFFT */
509  N2[t] = ceil(ths->sigma[t]*(ths->aN1[t]));
510 
511  /* N2 should be even */
512  if(N2[t]%2 != 0)
513  N2[t] = N2[t] +1;
514  }
515 
516  WINDOW_HELP_INIT
517 
518  if(ths->nnfft_flags & MALLOC_X)
519  ths->x = (double*)nfft_malloc(ths->d*ths->M_total*sizeof(double));
520  if(ths->nnfft_flags & MALLOC_F)
521  ths->f=(double _Complex*)nfft_malloc(ths->M_total*sizeof(double _Complex));
522 
523  if(ths->nnfft_flags & MALLOC_V)
524  ths->v = (double*)nfft_malloc(ths->d*ths->N_total*sizeof(double));
525  if(ths->nnfft_flags & MALLOC_F_HAT)
526  ths->f_hat = (double _Complex*)nfft_malloc(ths->N_total*sizeof(double _Complex));
527 
528  if(ths->nnfft_flags & PRE_LIN_PSI)
529  {
530  ths->K=(1U<< 10)*(ths->m+1);
531  ths->psi = (double*) nfft_malloc((ths->K+1)*ths->d*sizeof(double));
532  }
533 
534  if(ths->nnfft_flags & PRE_PSI)
535  ths->psi = (double*)nfft_malloc(ths->N_total*ths->d*(2*ths->m+2)*sizeof(double));
536 
537  if(ths->nnfft_flags & PRE_FULL_PSI)
538  {
539  for(t=0,lprod = 1; t<ths->d; t++)
540  lprod *= 2*ths->m+2;
541 
542  ths->psi = (double*)nfft_malloc(ths->N_total*lprod*sizeof(double));
543 
544  ths->psi_index_f = (int*) nfft_malloc(ths->N_total*sizeof(int));
545  ths->psi_index_g = (int*) nfft_malloc(ths->N_total*lprod*sizeof(int));
546  }
547 
548  ths->direct_plan = (nfft_plan*)nfft_malloc(sizeof(nfft_plan));
549 
550  nfft_init_guru(ths->direct_plan, ths->d, ths->aN1, ths->M_total, N2, m2,
551  nfft_flags, fftw_flags);
552 
553  ths->direct_plan->x = ths->x;
554  ths->direct_plan->f = ths->f;
555  ths->F = ths->direct_plan->f_hat;
556 
557  ths->mv_trafo = (void (*) (void* ))nnfft_trafo;
558  ths->mv_adjoint = (void (*) (void* ))nnfft_adjoint;
559 }
560 
561 void nnfft_init_guru(nnfft_plan *ths, int d, int N_total, int M_total, int *N, int *N1,
562  int m, unsigned nnfft_flags)
563 {
564  int t;
566  unsigned nfft_flags;
567  unsigned fftw_flags;
568 
569  ths->d= d;
570  ths->M_total= M_total;
571  ths->N_total= N_total;
572  ths->m= m;
573  ths->nnfft_flags= nnfft_flags;
574  fftw_flags= FFTW_ESTIMATE| FFTW_DESTROY_INPUT;
575  nfft_flags= PRE_PHI_HUT| MALLOC_F_HAT| FFTW_INIT| FFT_OUT_OF_PLACE;
576 
577  if(ths->nnfft_flags & PRE_PSI)
578  nfft_flags = nfft_flags | PRE_PSI;
579 
580  if(ths->nnfft_flags & PRE_FULL_PSI)
581  nfft_flags = nfft_flags | PRE_FULL_PSI;
582 
583  if(ths->nnfft_flags & PRE_LIN_PSI)
584  nfft_flags = nfft_flags | PRE_LIN_PSI;
585 
586  ths->N = (int*) nfft_malloc(ths->d*sizeof(int));
587  ths->N1 = (int*) nfft_malloc(ths->d*sizeof(int));
588 
589  for(t=0; t<d; t++) {
590  ths->N[t] = N[t];
591  ths->N1[t] = N1[t];
592  }
593  nnfft_init_help(ths,m,nfft_flags,fftw_flags);
594 }
595 
596 void nnfft_init(nnfft_plan *ths, int d, int N_total, int M_total, int *N)
597 {
598  int t;
600  unsigned nfft_flags;
601  unsigned fftw_flags;
602 
603  ths->d = d;
604  ths->M_total = M_total;
605  ths->N_total = N_total;
606 
607  /* m should be greater to get the same accuracy as the nfft */
608 /* Was soll dieser Ausdruck machen? Es handelt sich um eine Ganzzahl!
609 
610  WINDOW_HELP_ESTIMATE_m;
611 */
612 
613  ths->N = (int*) nfft_malloc(ths->d*sizeof(int));
614  ths->N1 = (int*) nfft_malloc(ths->d*sizeof(int));
615 
616  for(t=0; t<d; t++) {
617  ths->N[t] = N[t];
618 
619  /* the standard oversampling factor in the nnfft is 1.5 */
620  ths->N1[t] = ceil(1.5*ths->N[t]);
621 
622  /* N1 should be even */
623  if(ths->N1[t]%2 != 0)
624  ths->N1[t] = ths->N1[t] +1;
625  }
626 
627  ths->nnfft_flags=PRE_PSI| PRE_PHI_HUT| MALLOC_X| MALLOC_V| MALLOC_F_HAT| MALLOC_F;
628  nfft_flags= PRE_PSI| PRE_PHI_HUT| MALLOC_F_HAT| FFTW_INIT| FFT_OUT_OF_PLACE;
629 
630  fftw_flags= FFTW_ESTIMATE| FFTW_DESTROY_INPUT;
631 
632  nnfft_init_help(ths,ths->m,nfft_flags,fftw_flags);
633 }
634 
635 void nnfft_finalize(nnfft_plan *ths)
636 {
637 
638  nfft_finalize(ths->direct_plan);
639 
640  nfft_free(ths->direct_plan);
641 
642  nfft_free(ths->aN1);
643  nfft_free(ths->N);
644  nfft_free(ths->N1);
645 
646  if(ths->nnfft_flags & PRE_FULL_PSI)
647  {
648  nfft_free(ths->psi_index_g);
649  nfft_free(ths->psi_index_f);
650  nfft_free(ths->psi);
651  }
652 
653  if(ths->nnfft_flags & PRE_PSI)
654  nfft_free(ths->psi);
655 
656  if(ths->nnfft_flags & PRE_LIN_PSI)
657  nfft_free(ths->psi);
658 
659  if(ths->nnfft_flags & PRE_PHI_HUT)
660  nfft_free(ths->c_phi_inv);
661 
662  if(ths->nnfft_flags & MALLOC_F)
663  nfft_free(ths->f);
664 
665  if(ths->nnfft_flags & MALLOC_F_HAT)
666  nfft_free(ths->f_hat);
667 
668  if(ths->nnfft_flags & MALLOC_X)
669  nfft_free(ths->x);
670 
671  if(ths->nnfft_flags & MALLOC_V)
672  nfft_free(ths->v);
673 }

Generated on Tue Apr 30 2013 by Doxygen 1.8.1