tensor-0.1.0
 All Data Structures Namespaces Functions Variables Typedefs Enumerations Enumerator Groups Pages
fftw_common.cc
1 // -*- mode: c++; fill-column: 80; c-basic-offset: 2; indent-tabs-mode: nil -*-
2 /*
3  Copyright (c) 2013 Juan Jose Garcia Ripoll
4 
5  Tensor is free software; you can redistribute it and/or modify it
6  under the terms of the GNU Library General Public License as published
7  by the Free Software Foundation; either version 2 of the License, or
8  (at your option) any later version.
9 
10  This program is distributed in the hope that it will be useful,
11  but WITHOUT ANY WARRANTY; without even the implied warranty of
12  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  GNU Library General Public License for more details.
14 
15  You should have received a copy of the GNU General Public License along
16  with this program; if not, write to the Free Software Foundation, Inc.,
17  51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
18  */
19 
20 #include "fftw_common.hpp"
21 
22 
23 namespace tensor {
24 
25 // basic FFTW along all degrees of freedom
26 void
27 do_fftw(fftw_complex* pin, fftw_complex* pout, const Indices& dims, int direction) {
28  fftw_plan plan;
29  int d = dims.size();
30 
31  if (d == 1) {
32  plan = fftw_plan_dft_1d(dims[0], pin, pout, direction, FFTW_ESTIMATE);
33  } else {
34  int dimensions[d];
35  for (int i = 0; i < d; i++) {
36  dimensions[d - i - 1] = dims[i];
37  }
38  plan = fftw_plan_dft(d, dimensions, pin, pout, direction, FFTW_ESTIMATE);
39  }
40 
41  fftw_execute(plan);
42  fftw_destroy_plan(plan);
43 }
44 
45 void
46 do_fftw(fftw_complex* pin, fftw_complex* pout, index dim, const Indices& dims, int direction) {
47  // Collect the dimensions for the transform
48  fftw_iodim fft_dim;
49  fftw_iodim loop_dims[dims.size()-1];
50  index stride = 1;
51 
52  for (index i = 0; i < dim; i++) {
53  loop_dims[i].n = dims[i];
54  loop_dims[i].is = stride;
55  loop_dims[i].os = stride;
56  stride *= loop_dims[i].n;
57  }
58 
59  fft_dim.n = dims[dim];
60  fft_dim.is = stride;
61  fft_dim.os = stride;
62  stride *= fft_dim.n;
63 
64  for (index i = dim+1; i < dims.size(); i++) {
65  loop_dims[i-1].n = dims[i];
66  loop_dims[i-1].is = stride;
67  loop_dims[i-1].os = stride;
68  stride *= loop_dims[i-1].n;
69  }
70 
71  // Now make a Guru plan and execute it.
72  fftw_plan plan = fftw_plan_guru_dft(1, &fft_dim, dims.size()-1, loop_dims,
73  pin, pout, direction, FFTW_ESTIMATE);
74  fftw_execute(plan);
75  fftw_destroy_plan(plan);
76 }
77 
78 void
79 do_fftw(fftw_complex* pin, fftw_complex* pout, const Booleans& convert, const Indices& dims, int direction) {
80  // Collect the dimensions for the transform
81  fftw_iodim fft_dims[dims.size()];
82  fftw_iodim loop_dims[dims.size()];
83  index fft_size = 0;
84  index loop_size = 0;
85  index stride = 1;
86 
87  for (index i = 0; i < dims.size(); i++) {
88  if (convert[i] == false) {
89  loop_dims[loop_size].n = dims[i];
90  loop_dims[loop_size].is = stride;
91  loop_dims[loop_size].os = stride;
92  loop_size++;
93  } else {
94  fft_dims[fft_size].n = dims[i];
95  fft_dims[fft_size].is = stride;
96  fft_dims[fft_size].os = stride;
97  fft_size++;
98  }
99 
100  stride *= dims[i];
101  }
102 
103  // Now make a Guru plan and execute it.
104  fftw_plan plan = fftw_plan_guru_dft(fft_size, fft_dims, loop_size, loop_dims,
105  pin, pout, direction, FFTW_ESTIMATE);
106  fftw_execute(plan);
107  fftw_destroy_plan(plan);
108 }
109 
110 } // namespace tensor