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