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