NFFT 3.5.3alpha
construct_data_2d1d.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
37static void construct(char * file, int N, int M, int Z, fftw_complex *mem)
38{
39 int j,z; /* some variables */
40 double tmp; /* a placeholder */
41 nfft_plan my_plan; /* plan for the two dimensional nfft */
42 FILE* fp;
43
44 /* initialise my_plan */
45 nfft_init_2d(&my_plan,N,N,M/Z);
46
47 fp=fopen("knots.dat","r");
48
49 for(j=0;j<my_plan.M_total;j++)
50 {
51 fscanf(fp,"%le %le %le",&my_plan.x[2*j+0],&my_plan.x[2*j+1],&tmp);
52 }
53 fclose(fp);
54
55 fp=fopen(file,"w");
56
57 for(z=0;z<Z;z++) {
58 tmp = (double) z;
59
60 for(j=0;j<N*N;j++)
61 my_plan.f_hat[j] = mem[(z*N*N+N*N*Z/2+j)%(N*N*Z)];
62
63 if(my_plan.flags & PRE_PSI)
64 nfft_precompute_psi(&my_plan);
65
66 nfft_trafo(&my_plan);
67
68 for(j=0;j<my_plan.M_total;j++)
69 {
70 fprintf(fp,"%le %le %le %le %le\n",my_plan.x[2*j+0],my_plan.x[2*j+1],tmp/Z-0.5,
71 creal(my_plan.f[j]),cimag(my_plan.f[j]));
72 }
73 }
74 fclose(fp);
75
76 nfft_finalize(&my_plan);
77}
78
83static void fft(int N,int M,int Z, fftw_complex *mem)
84{
85 fftw_plan plan;
86 plan = fftw_plan_many_dft(1, &Z, N*N,
87 mem, NULL,
88 N*N, 1,
89 mem, NULL,
90 N*N,1 ,
91 FFTW_FORWARD, FFTW_ESTIMATE);
92
93 fftw_execute(plan); /* execute the fft */
94 fftw_destroy_plan(plan);
95}
96
101static void read_data(int N,int M,int Z, fftw_complex *mem)
102{
103 int i,z;
104 double real;
105 FILE* fin;
106 fin=fopen("input_f.dat","r");
107
108 for(z=0;z<Z;z++)
109 {
110 for(i=0;i<N*N;i++)
111 {
112 fscanf(fin,"%le ",&real );
113 mem[(z*N*N+N*N*Z/2+i)%(N*N*Z)]=real;
114 }
115 }
116 fclose(fin);
117}
118
119int main(int argc, char **argv)
120{
121 fftw_complex *mem;
122
123 if (argc <= 4) {
124 printf("usage: ./construct_data FILENAME N M Z\n");
125 return 1;
126 }
127
128 mem = (fftw_complex*) nfft_malloc(sizeof(fftw_complex) * atoi(argv[2]) * atoi(argv[2]) * atoi(argv[4]));
129
130 read_data(atoi(argv[2]),atoi(argv[3]),atoi(argv[4]), mem);
131
132 fft(atoi(argv[2]),atoi(argv[3]),atoi(argv[4]), mem);
133
134 construct(argv[1],atoi(argv[2]),atoi(argv[3]),atoi(argv[4]), mem);
135
136 nfft_free(mem);
137
138 return 1;
139}
140/* \} */
static void read_data(int N, int M, int Z, fftw_complex *mem)
read fills the memory with the file input_data_f.dat as the real part of f and with zeros for the ima...
static void construct(char *file, int N, int M, int Z, fftw_complex *mem)
construct makes an 2d-nfft for every slice
static void fft(int N, int M, int Z, fftw_complex *mem)
fft makes an 1D-ftt for every knot through all layers
#define PRE_PSI
Definition nfft3.h:191
Header file for the nfft3 library.
void * nfft_malloc(size_t n)
void nfft_free(void *p)