NFFT  3.4.1
construct_data_3d.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 
34 static void construct(char * file, int N, int M, int Z)
35 {
36  int j,k,l; /* some variables */
37  double real;
38  nfft_plan my_plan; /* plan for the three dimensional nfft */
39  FILE* fp,*fk;
40  int my_N[3],my_n[3]; /* to init the nfft */
41 
42 
43  /* initialise my_plan */
44  //nfft_init_3d(&my_plan,Z,N,N,M);
45  my_N[0]=Z; my_n[0]=ceil(Z*1.2);
46  my_N[1]=N; my_n[1]=ceil(N*1.2);
47  my_N[2]=N; my_n[2]=ceil(N*1.2);
48  nfft_init_guru(&my_plan, 3, my_N, M, my_n, 6,
51  FFTW_MEASURE| FFTW_DESTROY_INPUT);
52 
53  fp=fopen("knots.dat","r");
54 
55  for(j=0;j<M;j++)
56  fscanf(fp,"%le %le %le",&my_plan.x[3*(j)+1],
57  &my_plan.x[3*(j)+2],&my_plan.x[3*(j)+0]);
58 
59  fclose(fp);
60 
61  fp=fopen("input_f.dat","r");
62  fk=fopen(file,"w");
63 
64  for(l=0;l<Z;l++) {
65  for(j=0;j<N;j++)
66  {
67  for(k=0;k<N;k++)
68  {
69  //fscanf(fp,"%le ",&my_plan.f_hat[(N*N*(Z-l)+N*j+k+N*N*Z/2)%(N*N*Z)][0]);
70  fscanf(fp,"%le ",&real);
71  my_plan.f_hat[(N*N*l+N*j+k)] = real;
72  }
73  }
74  }
75 
76  if(my_plan.flags & PRE_PSI)
77  nfft_precompute_psi(&my_plan);
78 
79  nfft_trafo(&my_plan);
80 
81 
82  for(j=0;j<my_plan.M_total;j++)
83  fprintf(fk,"%le %le %le %le %le\n",my_plan.x[3*j+1],
84  my_plan.x[3*j+2],my_plan.x[3*j+0],creal(my_plan.f[j]),cimag(my_plan.f[j]));
85 
86 
87 
88  fclose(fk);
89  fclose(fp);
90 
91  nfft_finalize(&my_plan);
92 }
93 
94 int main(int argc, char **argv)
95 {
96  if (argc <= 4) {
97  printf("usage: ./construct_data FILENAME N M Z\n");
98  return 1;
99  }
100 
101  construct(argv[1], atoi(argv[2]),atoi(argv[3]),atoi(argv[4]));
102 
103  return 1;
104 }
105 /* \} */
fftw_complex * f_hat
Fourier coefficients.
Definition: nfft3.h:192
#define MALLOC_X
Definition: nfft3.h:199
#define MALLOC_F_HAT
Definition: nfft3.h:200
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
#define FFTW_INIT
Definition: nfft3.h:203
#define MALLOC_F
Definition: nfft3.h:201
NFFT_INT M_total
Total number of samples.
Definition: nfft3.h:192
#define FFT_OUT_OF_PLACE
Definition: nfft3.h:202
#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
#define PRE_PHI_HUT
Definition: nfft3.h:193
unsigned flags
Flags for precomputation, (de)allocation, and FFTW usage, default setting is PRE_PHI_HUT | PRE_PSI | ...
Definition: nfft3.h:192