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