NFFT  3.3.0
construct_data_2d.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 
21 #include "config.h"
22 
23 #include <stdlib.h>
24 #include <math.h>
25 #ifdef HAVE_COMPLEX_H
26 #include <complex.h>
27 #endif
28 
29 #include "nfft3.h"
30 
40 static void construct(char * file, int N, int M)
41 {
42  int j,k; /* some variables */
43  double real;
44  nfft_plan my_plan; /* plan for the two dimensional nfft */
45  FILE* fp;
46  FILE* fk;
47  FILE* fi;
48 
49  /* initialise my_plan */
50  nfft_init_2d(&my_plan,N,N,M);
51 
52  fp=fopen("knots.dat","r");
53 
54  for(j=0;j<my_plan.M_total;j++)
55  {
56  fscanf(fp,"%le %le ",&my_plan.x[2*j+0],&my_plan.x[2*j+1]);
57  }
58  fclose(fp);
59 
60  fi=fopen("input_f.dat","r");
61  fk=fopen(file,"w");
62 
63  for(j=0;j<N;j++)
64  {
65  for(k=0;k<N;k++) {
66  fscanf(fi,"%le ",&real);
67  my_plan.f_hat[(N*j+k)] = real;
68  }
69  }
70 
71  if(my_plan.flags & PRE_PSI)
72  nfft_precompute_psi(&my_plan);
73 
74  nfft_trafo(&my_plan);
75 
76  for(j=0;j<my_plan.M_total;j++)
77  {
78  fprintf(fk,"%le %le %le %le\n",my_plan.x[2*j+0],my_plan.x[2*j+1],creal(my_plan.f[j]),cimag(my_plan.f[j]));
79  }
80  fclose(fk);
81  fclose(fi);
82 
83  nfft_finalize(&my_plan);
84 }
85 
86 int main(int argc, char **argv)
87 {
88  if (argc <= 3) {
89  printf("usage: ./construct_data FILENAME N M\n");
90  return 1;
91  }
92 
93  construct(argv[1],atoi(argv[2]),atoi(argv[3]));
94 
95  return 1;
96 }
97 /* \} */
fftw_complex * f_hat
Fourier coefficients.
Definition: nfft3.h:194
static void construct(char *file, int N, int M)
construct makes an 2d-nfft
fftw_complex * f
Samples.
Definition: nfft3.h:194
data structure for an NFFT (nonequispaced fast Fourier transform) plan with double precision ...
Definition: nfft3.h:194
NFFT_INT M_total
Total number of samples.
Definition: nfft3.h:194
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