NFFT  3.3.0
construct_data_inh_2d1d.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 #include <limits.h>
25 #ifdef HAVE_COMPLEX_H
26 #include <complex.h>
27 #endif
28 
29 #include "nfft3.h"
30 
31 #ifndef MAX
32 #define MAX(a,b) (((a)>(b))?(a):(b))
33 #endif
34 
44 static void construct(char * file, int N, int M)
45 {
46  int j; /* some variables */
47  double real;
48  double w;
49  double time,min_time,max_time,min_inh,max_inh;
50  mri_inh_2d1d_plan my_plan;
51  FILE *fp,*fout,*fi,*finh,*ftime;
52  int my_N[3],my_n[3];
53  int flags = PRE_PHI_HUT| PRE_PSI |MALLOC_X| MALLOC_F_HAT|
54  MALLOC_F| FFTW_INIT| FFT_OUT_OF_PLACE|
55  FFTW_MEASURE| FFTW_DESTROY_INPUT;
56 
57  double Ts;
58  double W,T;
59  int N3;
60  int m=2;
61  double sigma = 1.25;
62 
63  ftime=fopen("readout_time.dat","r");
64  finh=fopen("inh.dat","r");
65 
66  min_time=INT_MAX; max_time=INT_MIN;
67  for(j=0;j<M;j++)
68  {
69  fscanf(ftime,"%le ",&time);
70  if(time<min_time)
71  min_time = time;
72  if(time>max_time)
73  max_time = time;
74  }
75 
76  fclose(ftime);
77 
78  Ts=(min_time+max_time)/2.0;
79 
80  min_inh=INT_MAX; max_inh=INT_MIN;
81  for(j=0;j<N*N;j++)
82  {
83  fscanf(finh,"%le ",&w);
84  if(w<min_inh)
85  min_inh = w;
86  if(w>max_inh)
87  max_inh = w;
88  }
89  fclose(finh);
90 
91 
92  N3=ceil((MAX(fabs(min_inh),fabs(max_inh))*(max_time-min_time)/2.0+m/(2*sigma))*4*sigma);
93  T=((max_time-min_time)/2.0)/(0.5-((double) m)/N3);
94  W=N3/T;
95 
96  my_N[0]=N; my_n[0]=ceil(N*sigma);
97  my_N[1]=N; my_n[1]=ceil(N*sigma);
98  my_N[2]=N3; my_n[2]=N3;
99 
100  /* initialise nfft */
101  mri_inh_2d1d_init_guru(&my_plan, my_N, M, my_n, m, sigma, flags,
102  FFTW_MEASURE| FFTW_DESTROY_INPUT);
103 
104  ftime=fopen("readout_time.dat","r");
105  fp=fopen("knots.dat","r");
106 
107  for(j=0;j<my_plan.M_total;j++)
108  {
109  fscanf(fp,"%le %le ",&my_plan.plan.x[2*j+0],&my_plan.plan.x[2*j+1]);
110  fscanf(ftime,"%le ",&my_plan.t[j]);
111  my_plan.t[j] = (my_plan.t[j]-Ts)/T;
112  }
113  fclose(fp);
114  fclose(ftime);
115 
116  finh=fopen("inh.dat","r");
117  for(j=0;j<N*N;j++)
118  {
119  fscanf(finh,"%le ",&my_plan.w[j]);
120  my_plan.w[j]/=W;
121  }
122  fclose(finh);
123 
124 
125  fi=fopen("input_f.dat","r");
126  for(j=0;j<N*N;j++)
127  {
128  fscanf(fi,"%le ",&real);
129  my_plan.f_hat[j] = real*cexp(2.0*_Complex_I*M_PI*Ts*my_plan.w[j]*W);
130  }
131 
132  if(my_plan.plan.flags & PRE_PSI)
133  nfft_precompute_psi(&my_plan.plan);
134 
135  mri_inh_2d1d_trafo(&my_plan);
136 
137  fout=fopen(file,"w");
138 
139  for(j=0;j<my_plan.M_total;j++)
140  {
141  fprintf(fout,"%le %le %le %le\n",my_plan.plan.x[2*j+0],my_plan.plan.x[2*j+1],creal(my_plan.f[j]),cimag(my_plan.f[j]));
142  }
143 
144  fclose(fout);
145 
146  mri_inh_2d1d_finalize(&my_plan);
147 }
148 
149 int main(int argc, char **argv)
150 {
151  if (argc <= 3) {
152  printf("usage: ./construct_data_inh_2d1d FILENAME N M\n");
153  return 1;
154  }
155 
156  construct(argv[1],atoi(argv[2]),atoi(argv[3]));
157 
158  return 1;
159 }
160 /* \} */
fftw_complex * f_hat
Fourier coefficients.
Definition: nfft3.h:527
fftw_complex * f
Samples.
Definition: nfft3.h:527
NFFT_INT M_total
Total number of samples.
Definition: nfft3.h:527
static void construct(char *file, int N, int M)
construct
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