NFFT Logo 3.2.3
mri3d/reconstruct_data_gridding.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: reconstruct_data_gridding.c 3896 2012-10-10 12:19:26Z tovo $ */
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 "nfft3util.h"
29 #include "nfft3.h"
30 
40 static void reconstruct(char* filename,int N,int M,int Z, int weight ,fftw_complex *mem)
41 {
42  int j,k,z; /* some variables */
43  double weights; /* store one weight temporary */
44  double tmp; /* tmp to read the obsolent z from the input file */
45  double real,imag; /* to read the real and imag part of a complex number */
46  nfft_plan my_plan; /* plan for the two dimensional nfft */
47  int my_N[2],my_n[2]; /* to init the nfft */
48  FILE* fin; /* input file */
49  FILE* fweight; /* input file for the weights */
50 
51  /* initialise my_plan */
52  my_N[0]=N; my_n[0]=ceil(N*1.2);
53  my_N[1]=N; my_n[1]=ceil(N*1.2);
54  nfft_init_guru(&my_plan, 2, my_N, M/Z, my_n, 6, PRE_PHI_HUT| PRE_PSI|
55  MALLOC_X| MALLOC_F_HAT| MALLOC_F|
56  FFTW_INIT| FFT_OUT_OF_PLACE,
57  FFTW_MEASURE| FFTW_DESTROY_INPUT);
58 
59  /* precompute lin psi if set */
60  if(my_plan.nfft_flags & PRE_LIN_PSI)
61  nfft_precompute_lin_psi(&my_plan);
62 
63  fin=fopen(filename,"r");
64 
65  for(z=0;z<Z;z++) {
66  fweight=fopen("weights.dat","r");
67  for(j=0;j<my_plan.M_total;j++)
68  {
69  fscanf(fweight,"%le ",&weights);
70  fscanf(fin,"%le %le %le %le %le",
71  &my_plan.x[2*j+0],&my_plan.x[2*j+1],&tmp,&real,&imag);
72  my_plan.f[j] = real + _Complex_I*imag;
73  if(weight)
74  my_plan.f[j] = my_plan.f[j] * weights;
75  }
76  fclose(fweight);
77 
78  /* precompute psi if set just one time because the knots equal each slice */
79  if(z==0 && my_plan.nfft_flags & PRE_PSI)
80  nfft_precompute_psi(&my_plan);
81 
82  /* precompute full psi if set just one time because the knots equal each slice */
83  if(z==0 && my_plan.nfft_flags & PRE_FULL_PSI)
84  nfft_precompute_full_psi(&my_plan);
85 
86  /* compute the adjoint nfft */
87  nfft_adjoint(&my_plan);
88 
89  for(k=0;k<my_plan.N_total;k++) {
90  /* write every slice in the memory.
91  here we make an fftshift direct */
92  mem[(Z*N*N/2+z*N*N+ k)%(Z*N*N)] = my_plan.f_hat[k];
93  }
94  }
95  fclose(fin);
96 
97  nfft_finalize(&my_plan);
98 }
99 
104 static void print(int N,int M,int Z, fftw_complex *mem)
105 {
106  int i,j;
107  FILE* fout_real;
108  FILE* fout_imag;
109  fout_real=fopen("output_real.dat","w");
110  fout_imag=fopen("output_imag.dat","w");
111 
112  for(i=0;i<Z;i++) {
113  for (j=0;j<N*N;j++) {
114  fprintf(fout_real,"%le ",creal(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
115  fprintf(fout_imag,"%le ",cimag(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
116  }
117  fprintf(fout_real,"\n");
118  fprintf(fout_imag,"\n");
119  }
120 
121  fclose(fout_real);
122  fclose(fout_imag);
123 }
124 
125 
126 int main(int argc, char **argv)
127 {
128  fftw_complex *mem;
129  fftw_plan plan;
130  int N,M,Z;
131 
132  if (argc <= 6) {
133  printf("usage: ./reconstruct_data_gridding FILENAME N M Z ITER WEIGHTS\n");
134  return 1;
135  }
136 
137  N=atoi(argv[2]);
138  M=atoi(argv[3]);
139  Z=atoi(argv[4]);
140 
141  /* Allocate memory to hold every slice in memory after the
142  2D-infft */
143  mem = (fftw_complex*) nfft_malloc(sizeof(fftw_complex) * atoi(argv[2]) * atoi(argv[2]) * atoi(argv[4]));
144 
145  /* Create plan for the 1d-ifft */
146  plan = fftw_plan_many_dft(1, &Z, N*N,
147  mem, NULL,
148  N*N, 1,
149  mem, NULL,
150  N*N,1 ,
151  FFTW_BACKWARD, FFTW_MEASURE);
152 
153  /* execute the 2d-nfft's */
154  reconstruct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]),atoi(argv[6]),mem);
155 
156  /* execute the 1d-fft's */
157  fftw_execute(plan);
158 
159  /* write the memory back in files */
160  print(N,M,Z, mem);
161 
162  /* free memory */
163  nfft_free(mem);
164 
165  return 1;
166 }
167 /* \} */

Generated on Tue Apr 30 2013 by Doxygen 1.8.1