NFFT  3.3.0
mri2d/reconstruct_data_gridding.c
1 /*
2  * Copyright (c) 2002, 2015 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$ */
20 #include "config.h"
21 
22 #include <stdlib.h>
23 #include <math.h>
24 #ifdef HAVE_COMPLEX_H
25 #include <complex.h>
26 #endif
27 
28 #include "nfft3.h"
29 
39 static void reconstruct(char* filename, int N, int M, int weight)
40 {
41  int j; /* some variables */
42  double weights; /* store one weight temporary */
43  double real,imag; /* to read the real and imag part of a complex number */
44  nfft_plan my_plan; /* plan for the two dimensional nfft */
45  FILE* fin; /* input file */
46  FILE* fweight; /* input file for the weights */
47  FILE *fout_real; /* output file */
48  FILE *fout_imag; /* output file */
49  int my_N[2],my_n[2];
50  int flags = PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT|
51  MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE|
52  FFTW_MEASURE| FFTW_DESTROY_INPUT;
53 
54  /* initialise nfft */
55  my_N[0]=N; my_n[0]=ceil(N*1.2);
56  my_N[1]=N; my_n[1]=ceil(N*1.2);
57  nfft_init_guru(&my_plan, 2, my_N, M, my_n, 6,flags,
58  FFTW_MEASURE| FFTW_DESTROY_INPUT);
59 
60  fin=fopen(filename,"r");
61 
62  fweight=fopen("weights.dat","r");
63  for(j=0;j<my_plan.M_total;j++)
64  {
65  fscanf(fweight,"%le ",&weights);
66  fscanf(fin,"%le %le %le %le",&my_plan.x[2*j+0],&my_plan.x[2*j+1],&real,&imag);
67  my_plan.f[j] = real + _Complex_I*imag;
68  if (weight)
69  my_plan.f[j] = my_plan.f[j] * weights;
70  }
71  fclose(fweight);
72 
73  /* precompute psi */
74  if(my_plan.flags & PRE_PSI)
75  nfft_precompute_psi(&my_plan);
76 
77  /* precompute full psi */
78  if(my_plan.flags & PRE_FULL_PSI)
79  nfft_precompute_full_psi(&my_plan);
80 
81 
82  /* compute the adjoint nfft */
83  nfft_adjoint(&my_plan);
84 
85  fout_real=fopen("output_real.dat","w");
86  fout_imag=fopen("output_imag.dat","w");
87 
88  for (j=0;j<N*N;j++) {
89  fprintf(fout_real,"%le ",creal(my_plan.f_hat[j]));
90  fprintf(fout_imag,"%le ",cimag(my_plan.f_hat[j]));
91  }
92 
93  fclose(fin);
94  fclose(fout_real);
95  fclose(fout_imag);
96 
97  nfft_finalize(&my_plan);
98 }
99 
100 
101 int main(int argc, char **argv)
102 {
103  if (argc <= 5) {
104  printf("usage: ./reconstruct_data_gridding FILENAME N M ITER WEIGHTS\n");
105  return 1;
106  }
107  reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[5]));
108 
109  return 1;
110 }
111 /* \} */
fftw_complex * f_hat
Fourier coefficients.
Definition: nfft3.h:194
fftw_complex * f
Samples.
Definition: nfft3.h:194
data structure for an NFFT (nonequispaced fast Fourier transform) plan with double precision ...
Definition: nfft3.h:194
NFFT_INT M_total
Total number of samples.
Definition: nfft3.h:194
static void reconstruct(char *filename, int N, int M, int weight)
reconstruct makes a 2d-adjoint-nfft
double * x
Nodes in time/spatial domain, size is doubles.
Definition: nfft3.h:194
unsigned flags
Flags for precomputation, (de)allocation, and FFTW usage, default setting is PRE_PHI_HUT | PRE_PSI | ...
Definition: nfft3.h:194