NFFT 3.5.3alpha
reconstruct_data_2d1d.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 iteration, int weight, fftw_complex *mem)
38{
39 int j,k,l,z; /* some variables */
40 double real,imag; /* to read the real and imag part of a complex number */
41 nfft_plan my_plan; /* plan for the two dimensional nfft */
42 solver_plan_complex my_iplan; /* plan for the two dimensional infft */
43 FILE* fin; /* input file */
44 int my_N[2],my_n[2]; /* to init the nfft */
45 double tmp, epsilon=0.0000003;/* tmp to read the obsolent z from the input file
46 epsilon is the break criterium for
47 the iteration */
48 unsigned infft_flags = CGNR | PRECOMPUTE_DAMP; /* flags for the infft */
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|
56 FFTW_MEASURE);
57
58 /* precompute lin psi if set */
59 if(my_plan.flags & PRE_LIN_PSI)
60 nfft_precompute_lin_psi(&my_plan);
61
62 /* set the flags for the infft*/
63 if (weight)
64 infft_flags = infft_flags | PRECOMPUTE_WEIGHT;
65
66 /* initialise my_iplan, advanced */
67 solver_init_advanced_complex(&my_iplan,(nfft_mv_plan_complex*)(&my_plan), infft_flags );
68
69 /* get the weights */
70 if(my_iplan.flags & PRECOMPUTE_WEIGHT)
71 {
72 fin=fopen("weights.dat","r");
73 for(j=0;j<my_plan.M_total;j++)
74 {
75 fscanf(fin,"%le ",&my_iplan.w[j]);
76 }
77 fclose(fin);
78 }
79
80 /* get the damping factors */
81 if(my_iplan.flags & PRECOMPUTE_DAMP)
82 {
83 for(j=0;j<N;j++){
84 for(k=0;k<N;k++) {
85 int j2= j-N/2;
86 int k2= k-N/2;
87 double r=sqrt(j2*j2+k2*k2);
88 if(r>(double) N/2)
89 my_iplan.w_hat[j*N+k]=0.0;
90 else
91 my_iplan.w_hat[j*N+k]=1.0;
92 }
93 }
94 }
95
96 /* open the input file */
97 fin=fopen(filename,"r");
98
99 /* For every Layer*/
100 for(z=0;z<Z;z++) {
101
102 /* read x,y,freal and fimag from the knots */
103 for(j=0;j<my_plan.M_total;j++)
104 {
105 fscanf(fin,"%le %le %le %le %le ",&my_plan.x[2*j+0],&my_plan.x[2*j+1], &tmp,
106 &real,&imag);
107 my_iplan.y[j] = real + _Complex_I*imag;
108 }
109
110 /* precompute psi if set just one time because the knots equal each plane */
111 if(z==0 && my_plan.flags & PRE_PSI)
112 nfft_precompute_psi(&my_plan);
113
114 /* precompute full psi if set just one time because the knots equal each plane */
115 if(z==0 && my_plan.flags & PRE_FULL_PSI)
116 nfft_precompute_full_psi(&my_plan);
117
118 /* init some guess */
119 for(k=0;k<my_plan.N_total;k++)
120 my_iplan.f_hat_iter[k]=0.0;
121
122 /* inverse trafo */
123 solver_before_loop_complex(&my_iplan);
124 for(l=0;l<iteration;l++)
125 {
126 /* break if dot_r_iter is smaller than epsilon*/
127 if(my_iplan.dot_r_iter<epsilon)
128 break;
129 fprintf(stderr,"%e, %i of %i\n",sqrt(my_iplan.dot_r_iter),
130 iteration*z+l+1,iteration*Z);
131 solver_loop_one_step_complex(&my_iplan);
132 }
133 for(k=0;k<my_plan.N_total;k++) {
134 /* write every slice in the memory.
135 here we make an fftshift direct */
136 mem[(Z*N*N/2+z*N*N+ k)%(Z*N*N)] = my_iplan.f_hat_iter[k];
137 }
138 }
139
140 fclose(fin);
141
142 /* finalize the infft */
143 solver_finalize_complex(&my_iplan);
144
145 /* finalize the nfft */
146 nfft_finalize(&my_plan);
147}
148
153static void print(int N,int M,int Z, fftw_complex *mem)
154{
155 int i,j;
156 FILE* fout_real;
157 FILE* fout_imag;
158 fout_real=fopen("output_real.dat","w");
159 fout_imag=fopen("output_imag.dat","w");
160
161 for(i=0;i<Z;i++) {
162 for (j=0;j<N*N;j++) {
163 fprintf(fout_real,"%le ",creal(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
164 fprintf(fout_imag,"%le ",cimag(mem[(Z*N*N/2+i*N*N+ j)%(Z*N*N)]) /Z);
165 }
166 fprintf(fout_real,"\n");
167 fprintf(fout_imag,"\n");
168 }
169
170 fclose(fout_real);
171 fclose(fout_imag);
172}
173
174int main(int argc, char **argv)
175{
176 fftw_complex *mem;
177 fftw_plan plan;
178 int N,M,Z;
179
180 if (argc <= 6) {
181 printf("usage: ./reconstruct FILENAME N M Z ITER WEIGHTS\n");
182 return 1;
183 }
184
185 N=atoi(argv[2]);
186 M=atoi(argv[3]);
187 Z=atoi(argv[4]);
188
189 /* Allocate memory to hold every layer in memory after the
190 2D-infft */
191 mem = (fftw_complex*) nfft_malloc(sizeof(fftw_complex) * atoi(argv[2]) * atoi(argv[2]) * atoi(argv[4]));
192
193 /* Create plan for the 1d-ifft */
194 plan = fftw_plan_many_dft(1, &Z, N*N,
195 mem, NULL,
196 N*N, 1,
197 mem, NULL,
198 N*N,1 ,
199 FFTW_BACKWARD, FFTW_MEASURE);
200
201 /* execute the 2d-infft's */
202 reconstruct(argv[1],N,M,Z,atoi(argv[5]),atoi(argv[6]),mem);
203
204 /* execute the 1d-fft's */
205 fftw_execute(plan);
206
207 /* write the memory back in files */
208 print(N,M,Z, mem);
209
210 /* free memory */
211 nfft_free(mem);
212 fftw_destroy_plan(plan);
213 return 1;
214}
215/* \} */
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....
static void reconstruct(char *filename, int N, int M, int Z, int iteration, int weight, fftw_complex *mem)
reconstruct makes an inverse 2d-nfft for every slice
#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
#define CGNR
Definition nfft3.h:808
#define PRECOMPUTE_DAMP
Definition nfft3.h:812
#define PRECOMPUTE_WEIGHT
Definition nfft3.h:811
Header file for the nfft3 library.
void * nfft_malloc(size_t n)
void nfft_free(void *p)
data structure for an inverse NFFT plan with double precision
Definition nfft3.h:802
double * w
weighting factors
Definition nfft3.h:802
unsigned flags
iteration type
Definition nfft3.h:802
double * w_hat
damping factors
Definition nfft3.h:802
double dot_r_iter
weighted dotproduct of r_iter
Definition nfft3.h:802
fftw_complex * y
right hand side, samples
Definition nfft3.h:802
fftw_complex * f_hat_iter
iterative solution
Definition nfft3.h:802