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