NFFT Logo 3.2.3
construct_data_2d1d.c
1 /*
2  * Copyright (c) 2002, 2012 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: construct_data_2d1d.c 3896 2012-10-10 12:19:26Z tovo $ */
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 "nfft3util.h"
29 #include "nfft3.h"
30 
40 static void construct(char * file, int N, int M, int Z, fftw_complex *mem)
41 {
42  int j,z; /* some variables */
43  double tmp; /* a placeholder */
44  nfft_plan my_plan; /* plan for the two dimensional nfft */
45  FILE* fp;
46 
47  /* initialise my_plan */
48  nfft_init_2d(&my_plan,N,N,M/Z);
49 
50  fp=fopen("knots.dat","r");
51 
52  for(j=0;j<my_plan.M_total;j++)
53  {
54  fscanf(fp,"%le %le %le",&my_plan.x[2*j+0],&my_plan.x[2*j+1],&tmp);
55  }
56  fclose(fp);
57 
58  fp=fopen(file,"w");
59 
60  for(z=0;z<Z;z++) {
61  tmp = (double) z;
62 
63  for(j=0;j<N*N;j++)
64  my_plan.f_hat[j] = mem[(z*N*N+N*N*Z/2+j)%(N*N*Z)];
65 
66  if(my_plan.nfft_flags & PRE_PSI)
67  nfft_precompute_psi(&my_plan);
68 
69  nfft_trafo(&my_plan);
70 
71  for(j=0;j<my_plan.M_total;j++)
72  {
73  fprintf(fp,"%le %le %le %le %le\n",my_plan.x[2*j+0],my_plan.x[2*j+1],tmp/Z-0.5,
74  creal(my_plan.f[j]),cimag(my_plan.f[j]));
75  }
76  }
77  fclose(fp);
78 
79  nfft_finalize(&my_plan);
80 }
81 
86 static void fft(int N,int M,int Z, fftw_complex *mem)
87 {
88  fftw_plan plan;
89  plan = fftw_plan_many_dft(1, &Z, N*N,
90  mem, NULL,
91  N*N, 1,
92  mem, NULL,
93  N*N,1 ,
94  FFTW_FORWARD, FFTW_ESTIMATE);
95 
96  fftw_execute(plan); /* execute the fft */
97  fftw_destroy_plan(plan);
98 }
99 
104 static void read_data(int N,int M,int Z, fftw_complex *mem)
105 {
106  int i,z;
107  double real;
108  FILE* fin;
109  fin=fopen("input_f.dat","r");
110 
111  for(z=0;z<Z;z++)
112  {
113  for(i=0;i<N*N;i++)
114  {
115  fscanf(fin,"%le ",&real );
116  mem[(z*N*N+N*N*Z/2+i)%(N*N*Z)]=real;
117  }
118  }
119  fclose(fin);
120 }
121 
122 int main(int argc, char **argv)
123 {
124  fftw_complex *mem;
125 
126  if (argc <= 4) {
127  printf("usage: ./construct_data FILENAME N M Z\n");
128  return 1;
129  }
130 
131  mem = (fftw_complex*) nfft_malloc(sizeof(fftw_complex) * atoi(argv[2]) * atoi(argv[2]) * atoi(argv[4]));
132 
133  read_data(atoi(argv[2]),atoi(argv[3]),atoi(argv[4]), mem);
134 
135  fft(atoi(argv[2]),atoi(argv[3]),atoi(argv[4]), mem);
136 
137  construct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]), mem);
138 
139  nfft_free(mem);
140 
141  return 1;
142 }
143 /* \} */

Generated on Tue Apr 30 2013 by Doxygen 1.8.1