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