Mercurial > octave
annotate liboctave/numeric/oct-fftw.cc @ 20955:77f5591878bf
maint: Use '! expr' rather than '!expr' to conform to coding guidelines.
* dialog.h, documentation-dock-widget.cc, files-dock-widget.cc,
find-files-dialog.cc, file-editor-tab.cc, file-editor.cc, find-dialog.cc,
octave-qscintilla.cc, main-window.cc, parser.cc, resource-manager.cc,
workspace-view.cc, data.cc, dlmread.cc, gl-render.cc, gl2ps-renderer.cc,
graphics.cc, graphics.in.h, ls-hdf5.cc, ls-mat5.cc, ls-oct-binary.cc, lsode.cc,
mappers.cc, pt-jit.cc, regexp.cc, spparms.cc, symtab.h, utils.cc, zfstream.cc,
__eigs__.cc, __glpk__.cc, __init_fltk__.cc, ccolamd.cc, colamd.cc,
ov-base-diag.cc, ov-base-int.cc, ov-base-sparse.cc, ov-bool-mat.cc,
ov-bool-sparse.cc, ov-bool.cc, ov-class.cc, ov-cx-sparse.cc, ov-fcn-handle.cc,
ov-fcn-inline.cc, ov-java.cc, ov-perm.cc, ov-re-sparse.cc, ov-str-mat.cc,
ov-struct.cc, ov.cc, pt-mat.cc, Array.cc, Array.h, CMatrix.cc, CSparse.cc,
MatrixType.cc, PermMatrix.cc, Sparse.h, dMatrix.cc, dSparse.cc, fCMatrix.cc,
fMatrix.cc, idx-vector.cc, CollocWt.h, SparseCmplxLU.cc, SparseCmplxQR.cc,
SparseQR.cc, SparsedbleLU.cc, base-qr.cc, eigs-base.cc, oct-fftw.cc,
randmtzig.c, sparse-dmsolve.cc, kpse.cc, lo-regexp.cc, oct-locbuf.h,
url-transfer.cc, url-transfer.h, bitset.m, interp2.m, __isequal__.m,
inpolygon.m, questdlg.m, help.m, compare_versions.m, substruct.m,
configure_make.m, whitebg.m, __marching_cube__.m, struct2hdl.m, polyfit.m,
spline.m, unique.m, svds.m, ellipke.m, binoinv.m, hygepdf.m, nbininv.m,
poissinv.m, tcdf.m, unidcdf.m, unidpdf.m, dec2base.m, assert.m, weekday.m,
mkoctfile.in.cc:
maint: Use '! expr' rather than '!expr' to conform to coding guidelines.
author | Rik <rik@octave.org> |
---|---|
date | Sun, 20 Dec 2015 10:15:02 -0800 |
parents | 4197fc428c7d |
children | 7cac4e7458f2 |
rev | line source |
---|---|
3828 | 1 /* |
2 | |
19697
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
17769
diff
changeset
|
3 Copyright (C) 2001-2015 John W. Eaton |
7017 | 4 |
3828 | 5 This file is part of Octave. |
6 | |
7 Octave is free software; you can redistribute it and/or modify it | |
8 under the terms of the GNU General Public License as published by the | |
7016 | 9 Free Software Foundation; either version 3 of the License, or (at your |
10 option) any later version. | |
3828 | 11 |
12 Octave is distributed in the hope that it will be useful, but WITHOUT | |
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License | |
15 for more details. | |
16 | |
17 You should have received a copy of the GNU General Public License | |
7016 | 18 along with Octave; see the file COPYING. If not, see |
19 <http://www.gnu.org/licenses/>. | |
3828 | 20 |
21 */ | |
22 | |
23 #ifdef HAVE_CONFIG_H | |
24 #include <config.h> | |
25 #endif | |
26 | |
9523
0ce82753dd72
more configure changes for libraries
John W. Eaton <jwe@octave.org>
parents:
9516
diff
changeset
|
27 #if defined (HAVE_FFTW) |
3828 | 28 |
4775 | 29 #include <iostream> |
30 #include <vector> | |
31 | |
4786 | 32 #include "lo-error.h" |
3828 | 33 #include "oct-fftw.h" |
4786 | 34 #include "quit.h" |
8377
25bc2d31e1bf
improve OCTAVE_LOCAL_BUFFER
Jaroslav Hajek <highegg@gmail.com>
parents:
7789
diff
changeset
|
35 #include "oct-locbuf.h" |
13983
7dd7cccf0757
clean up memory allocated for singletons before exit
John W. Eaton <jwe@octave.org>
parents:
11586
diff
changeset
|
36 #include "singleton-cleanup.h" |
3828 | 37 |
15960
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
38 #if defined (HAVE_FFTW3_THREADS) |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
39 #include "nproc.h" |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
40 #endif |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
41 |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
42 octave_fftw_planner *octave_fftw_planner::instance = 0; |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
43 |
14078
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
44 // Helper class to create and cache FFTW plans for both 1D and |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
45 // 2D. This implementation defaults to using FFTW_ESTIMATE to create |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
46 // the plans, which in theory is suboptimal, but provides quite |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
47 // reasonable performance in practice. |
4773 | 48 |
14078
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
49 // Also note that if FFTW_ESTIMATE is not used then the planner in FFTW3 |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
50 // will destroy the input and output arrays. We must, therefore, create a |
6228 | 51 // temporary input array with the same size and 16-byte alignment as |
14078
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
52 // the original array when using a different planner strategy. |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
53 // Note that we also use any wisdom that is available, either in a |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
54 // FFTW3 system wide file or as supplied by the user. |
4773 | 55 |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
56 // FIXME: if we can ensure 16 byte alignment in Array<T> |
4809 | 57 // (<T> *data) the FFTW3 can use SIMD instructions for further |
58 // acceleration. | |
4773 | 59 |
14078
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
60 // Note that it is profitable to store the FFTW3 plans, for small FFTs. |
3828 | 61 |
4808 | 62 octave_fftw_planner::octave_fftw_planner (void) |
11509
fc35194006d6
oct-fftw.cc: more constructor tweaks
John W. Eaton <jwe@octave.org>
parents:
11501
diff
changeset
|
63 : meth (ESTIMATE), rplan (0), rd (0), rs (0), rr (0), rh (0), rn (), |
11501
331fcc41ca23
data member initialization fixes
John W. Eaton <jwe@octave.org>
parents:
10350
diff
changeset
|
64 rsimd_align (false) |
3828 | 65 { |
66 plan[0] = plan[1] = 0; | |
4773 | 67 d[0] = d[1] = s[0] = s[1] = r[0] = r[1] = h[0] = h[1] = 0; |
4808 | 68 simd_align[0] = simd_align[1] = false; |
5044 | 69 inplace[0] = inplace[1] = false; |
4808 | 70 n[0] = n[1] = dim_vector (); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
71 |
15960
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
72 #if defined (HAVE_FFTW3_THREADS) |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
73 int init_ret = fftw_init_threads (); |
20955
77f5591878bf
maint: Use '! expr' rather than '!expr' to conform to coding guidelines.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
74 if (! init_ret) |
15960
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
75 (*current_liboctave_error_handler) ("Error initializing FFTW threads"); |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
76 //Use number of processors available to the current process |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
77 //This can be later changed with fftw ("threads", nthreads) |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
78 nthreads = num_processors (NPROC_CURRENT); |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
79 fftw_plan_with_nthreads (nthreads); |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
80 #endif |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
81 |
4809 | 82 // If we have a system wide wisdom file, import it. |
4808 | 83 fftw_import_system_wisdom (); |
3828 | 84 } |
85 | |
14078
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
86 octave_fftw_planner::~octave_fftw_planner (void) |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
87 { |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
88 fftw_plan *plan_p; |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
89 |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
90 plan_p = &rplan; |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
91 if (*plan_p) |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
92 fftw_destroy_plan (*plan_p); |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
93 |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
94 plan_p = &plan[0]; |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
95 if (*plan_p) |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
96 fftw_destroy_plan (*plan_p); |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
97 |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
98 plan_p = &plan[1]; |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
99 if (*plan_p) |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
100 fftw_destroy_plan (*plan_p); |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
101 } |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
102 |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
103 bool |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
104 octave_fftw_planner::instance_ok (void) |
6228 | 105 { |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
106 bool retval = true; |
6228 | 107 |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
108 if (! instance) |
13983
7dd7cccf0757
clean up memory allocated for singletons before exit
John W. Eaton <jwe@octave.org>
parents:
11586
diff
changeset
|
109 { |
7dd7cccf0757
clean up memory allocated for singletons before exit
John W. Eaton <jwe@octave.org>
parents:
11586
diff
changeset
|
110 instance = new octave_fftw_planner (); |
7dd7cccf0757
clean up memory allocated for singletons before exit
John W. Eaton <jwe@octave.org>
parents:
11586
diff
changeset
|
111 |
7dd7cccf0757
clean up memory allocated for singletons before exit
John W. Eaton <jwe@octave.org>
parents:
11586
diff
changeset
|
112 if (instance) |
7dd7cccf0757
clean up memory allocated for singletons before exit
John W. Eaton <jwe@octave.org>
parents:
11586
diff
changeset
|
113 singleton_cleanup_list::add (cleanup_instance); |
7dd7cccf0757
clean up memory allocated for singletons before exit
John W. Eaton <jwe@octave.org>
parents:
11586
diff
changeset
|
114 } |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
115 |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
116 if (! instance) |
6228 | 117 { |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
118 (*current_liboctave_error_handler) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
119 ("unable to create octave_fftw_planner object!"); |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
120 |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
121 retval = false; |
6228 | 122 } |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
123 |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
124 return retval; |
6228 | 125 } |
126 | |
4808 | 127 #define CHECK_SIMD_ALIGNMENT(x) \ |
6482 | 128 (((reinterpret_cast<ptrdiff_t> (x)) & 0xF) == 0) |
4808 | 129 |
3828 | 130 fftw_plan |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
131 octave_fftw_planner::do_create_plan (int dir, const int rank, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
132 const dim_vector dims, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
133 octave_idx_type howmany, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
134 octave_idx_type stride, |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
135 octave_idx_type dist, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
136 const Complex *in, Complex *out) |
3828 | 137 { |
4773 | 138 int which = (dir == FFTW_FORWARD) ? 0 : 1; |
3828 | 139 fftw_plan *cur_plan_p = &plan[which]; |
140 bool create_new_plan = false; | |
4808 | 141 bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out); |
5044 | 142 bool ioinplace = (in == out); |
3828 | 143 |
4809 | 144 // Don't create a new plan if we have a non SIMD plan already but |
145 // can do SIMD. This prevents endlessly recreating plans if we | |
146 // change the alignment. | |
147 | |
4783 | 148 if (plan[which] == 0 || d[which] != dist || s[which] != stride |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
149 || r[which] != rank || h[which] != howmany |
5044 | 150 || ioinplace != inplace[which] |
20955
77f5591878bf
maint: Use '! expr' rather than '!expr' to conform to coding guidelines.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
151 || ((ioalign != simd_align[which]) ? ! ioalign : false)) |
4773 | 152 create_new_plan = true; |
153 else | |
4809 | 154 { |
155 // We still might not have the same shape of array. | |
156 | |
157 for (int i = 0; i < rank; i++) | |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
158 if (dims(i) != n[which](i)) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
159 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
160 create_new_plan = true; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
161 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
162 } |
4809 | 163 } |
3828 | 164 |
165 if (create_new_plan) | |
166 { | |
4773 | 167 d[which] = dist; |
168 s[which] = stride; | |
169 r[which] = rank; | |
170 h[which] = howmany; | |
4808 | 171 simd_align[which] = ioalign; |
5044 | 172 inplace[which] = ioinplace; |
4773 | 173 n[which] = dims; |
174 | |
6228 | 175 // Note reversal of dimensions for column major storage in FFTW. |
176 octave_idx_type nn = 1; | |
177 OCTAVE_LOCAL_BUFFER (int, tmp, rank); | |
178 | |
179 for (int i = 0, j = rank-1; i < rank; i++, j--) | |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
180 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
181 tmp[i] = dims(j); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
182 nn *= dims(j); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
183 } |
6228 | 184 |
185 int plan_flags = 0; | |
186 bool plan_destroys_in = true; | |
187 | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
188 switch (meth) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
189 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
190 case UNKNOWN: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
191 case ESTIMATE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
192 plan_flags |= FFTW_ESTIMATE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
193 plan_destroys_in = false; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
194 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
195 case MEASURE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
196 plan_flags |= FFTW_MEASURE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
197 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
198 case PATIENT: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
199 plan_flags |= FFTW_PATIENT; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
200 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
201 case EXHAUSTIVE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
202 plan_flags |= FFTW_EXHAUSTIVE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
203 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
204 case HYBRID: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
205 if (nn < 8193) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
206 plan_flags |= FFTW_MEASURE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
207 else |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
208 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
209 plan_flags |= FFTW_ESTIMATE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
210 plan_destroys_in = false; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
211 } |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
212 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
213 } |
6228 | 214 |
4808 | 215 if (ioalign) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
216 plan_flags &= ~FFTW_UNALIGNED; |
4808 | 217 else |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
218 plan_flags |= FFTW_UNALIGNED; |
4808 | 219 |
3828 | 220 if (*cur_plan_p) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
221 fftw_destroy_plan (*cur_plan_p); |
3828 | 222 |
6228 | 223 if (plan_destroys_in) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
224 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
225 // Create matrix with the same size and 16-byte alignment as input |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
226 OCTAVE_LOCAL_BUFFER (Complex, itmp, nn * howmany + 32); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
227 itmp = reinterpret_cast<Complex *> |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
228 (((reinterpret_cast<ptrdiff_t>(itmp) + 15) & ~ 0xF) + |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
229 ((reinterpret_cast<ptrdiff_t> (in)) & 0xF)); |
4809 | 230 |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
231 *cur_plan_p = |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
232 fftw_plan_many_dft (rank, tmp, howmany, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
233 reinterpret_cast<fftw_complex *> (itmp), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
234 0, stride, dist, |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
235 reinterpret_cast<fftw_complex *> (out), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
236 0, stride, dist, dir, plan_flags); |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
237 } |
6228 | 238 else |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
239 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
240 *cur_plan_p = |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
241 fftw_plan_many_dft (rank, tmp, howmany, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
242 reinterpret_cast<fftw_complex *> (const_cast<Complex *> (in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
243 0, stride, dist, |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
244 reinterpret_cast<fftw_complex *> (out), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
245 0, stride, dist, dir, plan_flags); |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
246 } |
3828 | 247 |
248 if (*cur_plan_p == 0) | |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
249 (*current_liboctave_error_handler) ("Error creating fftw plan"); |
3828 | 250 } |
251 | |
252 return *cur_plan_p; | |
253 } | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
254 |
4773 | 255 fftw_plan |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
256 octave_fftw_planner::do_create_plan (const int rank, const dim_vector dims, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
257 octave_idx_type howmany, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
258 octave_idx_type stride, |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
259 octave_idx_type dist, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
260 const double *in, Complex *out) |
3828 | 261 { |
4773 | 262 fftw_plan *cur_plan_p = &rplan; |
3828 | 263 bool create_new_plan = false; |
4808 | 264 bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out); |
3828 | 265 |
4809 | 266 // Don't create a new plan if we have a non SIMD plan already but |
267 // can do SIMD. This prevents endlessly recreating plans if we | |
268 // change the alignment. | |
269 | |
4783 | 270 if (rplan == 0 || rd != dist || rs != stride || rr != rank |
20955
77f5591878bf
maint: Use '! expr' rather than '!expr' to conform to coding guidelines.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
271 || rh != howmany || ((ioalign != rsimd_align) ? ! ioalign : false)) |
4773 | 272 create_new_plan = true; |
273 else | |
4809 | 274 { |
275 // We still might not have the same shape of array. | |
276 | |
277 for (int i = 0; i < rank; i++) | |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
278 if (dims(i) != rn(i)) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
279 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
280 create_new_plan = true; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
281 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
282 } |
4809 | 283 } |
3828 | 284 |
285 if (create_new_plan) | |
286 { | |
4773 | 287 rd = dist; |
288 rs = stride; | |
289 rr = rank; | |
290 rh = howmany; | |
4808 | 291 rsimd_align = ioalign; |
4773 | 292 rn = dims; |
293 | |
6228 | 294 // Note reversal of dimensions for column major storage in FFTW. |
295 octave_idx_type nn = 1; | |
296 OCTAVE_LOCAL_BUFFER (int, tmp, rank); | |
297 | |
298 for (int i = 0, j = rank-1; i < rank; i++, j--) | |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
299 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
300 tmp[i] = dims(j); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
301 nn *= dims(j); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
302 } |
6228 | 303 |
304 int plan_flags = 0; | |
305 bool plan_destroys_in = true; | |
306 | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
307 switch (meth) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
308 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
309 case UNKNOWN: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
310 case ESTIMATE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
311 plan_flags |= FFTW_ESTIMATE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
312 plan_destroys_in = false; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
313 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
314 case MEASURE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
315 plan_flags |= FFTW_MEASURE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
316 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
317 case PATIENT: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
318 plan_flags |= FFTW_PATIENT; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
319 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
320 case EXHAUSTIVE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
321 plan_flags |= FFTW_EXHAUSTIVE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
322 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
323 case HYBRID: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
324 if (nn < 8193) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
325 plan_flags |= FFTW_MEASURE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
326 else |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
327 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
328 plan_flags |= FFTW_ESTIMATE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
329 plan_destroys_in = false; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
330 } |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
331 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
332 } |
6228 | 333 |
4808 | 334 if (ioalign) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
335 plan_flags &= ~FFTW_UNALIGNED; |
4808 | 336 else |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
337 plan_flags |= FFTW_UNALIGNED; |
4808 | 338 |
3828 | 339 if (*cur_plan_p) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
340 fftw_destroy_plan (*cur_plan_p); |
3828 | 341 |
6228 | 342 if (plan_destroys_in) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
343 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
344 // Create matrix with the same size and 16-byte alignment as input |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
345 OCTAVE_LOCAL_BUFFER (double, itmp, nn + 32); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
346 itmp = reinterpret_cast<double *> |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
347 (((reinterpret_cast<ptrdiff_t>(itmp) + 15) & ~ 0xF) + |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
348 ((reinterpret_cast<ptrdiff_t> (in)) & 0xF)); |
4809 | 349 |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
350 *cur_plan_p = |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
351 fftw_plan_many_dft_r2c (rank, tmp, howmany, itmp, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
352 0, stride, dist, |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
353 reinterpret_cast<fftw_complex *> (out), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
354 0, stride, dist, plan_flags); |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
355 } |
6228 | 356 else |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
357 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
358 *cur_plan_p = |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
359 fftw_plan_many_dft_r2c (rank, tmp, howmany, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
360 (const_cast<double *> (in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
361 0, stride, dist, |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
362 reinterpret_cast<fftw_complex *> (out), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
363 0, stride, dist, plan_flags); |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
364 } |
3828 | 365 |
366 if (*cur_plan_p == 0) | |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
367 (*current_liboctave_error_handler) ("Error creating fftw plan"); |
3828 | 368 } |
369 | |
370 return *cur_plan_p; | |
371 } | |
372 | |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
373 octave_fftw_planner::FftwMethod |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
374 octave_fftw_planner::do_method (void) |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
375 { |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
376 return meth; |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
377 } |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
378 |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
379 octave_fftw_planner::FftwMethod |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
380 octave_fftw_planner::do_method (FftwMethod _meth) |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
381 { |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
382 FftwMethod ret = meth; |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
383 if (_meth == ESTIMATE || _meth == MEASURE |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
384 || _meth == PATIENT || _meth == EXHAUSTIVE |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
385 || _meth == HYBRID) |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
386 { |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
387 if (meth != _meth) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
388 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
389 meth = _meth; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
390 if (rplan) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
391 fftw_destroy_plan (rplan); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
392 if (plan[0]) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
393 fftw_destroy_plan (plan[0]); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
394 if (plan[1]) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
395 fftw_destroy_plan (plan[1]); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
396 rplan = plan[0] = plan[1] = 0; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
397 } |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
398 } |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
399 else |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
400 ret = UNKNOWN; |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
401 return ret; |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
402 } |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
403 |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
404 octave_float_fftw_planner *octave_float_fftw_planner::instance = 0; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
405 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
406 octave_float_fftw_planner::octave_float_fftw_planner (void) |
11509
fc35194006d6
oct-fftw.cc: more constructor tweaks
John W. Eaton <jwe@octave.org>
parents:
11501
diff
changeset
|
407 : meth (ESTIMATE), rplan (0), rd (0), rs (0), rr (0), rh (0), rn (), |
11501
331fcc41ca23
data member initialization fixes
John W. Eaton <jwe@octave.org>
parents:
10350
diff
changeset
|
408 rsimd_align (false) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
409 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
410 plan[0] = plan[1] = 0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
411 d[0] = d[1] = s[0] = s[1] = r[0] = r[1] = h[0] = h[1] = 0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
412 simd_align[0] = simd_align[1] = false; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
413 inplace[0] = inplace[1] = false; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
414 n[0] = n[1] = dim_vector (); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
415 |
15960
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
416 #if defined (HAVE_FFTW3F_THREADS) |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
417 int init_ret = fftwf_init_threads (); |
20955
77f5591878bf
maint: Use '! expr' rather than '!expr' to conform to coding guidelines.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
418 if (! init_ret) |
15960
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
419 (*current_liboctave_error_handler) ("Error initializing FFTW3F threads"); |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
420 //Use number of processors available to the current process |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
421 //This can be later changed with fftw ("threads", nthreads) |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
422 nthreads = num_processors (NPROC_CURRENT); |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
423 fftwf_plan_with_nthreads (nthreads); |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
424 #endif |
bde7731b2b83
added FFTW multithreaded library support
Andreas Weber <andy.weber.aw@gmail.com>
parents:
15271
diff
changeset
|
425 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
426 // If we have a system wide wisdom file, import it. |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
427 fftwf_import_system_wisdom (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
428 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
429 |
14078
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
430 octave_float_fftw_planner::~octave_float_fftw_planner (void) |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
431 { |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
432 fftwf_plan *plan_p; |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
433 |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
434 plan_p = &rplan; |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
435 if (*plan_p) |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
436 fftwf_destroy_plan (*plan_p); |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
437 |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
438 plan_p = &plan[0]; |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
439 if (*plan_p) |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
440 fftwf_destroy_plan (*plan_p); |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
441 |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
442 plan_p = &plan[1]; |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
443 if (*plan_p) |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
444 fftwf_destroy_plan (*plan_p); |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
445 } |
941d19370065
Cleanup FFTW wisdom plans in class destructor and prevent a memory leak.
Rik <octave@nomad.inbox5.com>
parents:
13994
diff
changeset
|
446 |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
447 bool |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
448 octave_float_fftw_planner::instance_ok (void) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
449 { |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
450 bool retval = true; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
451 |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
452 if (! instance) |
13994
99f039289e95
also clean up float fftw planner
John W. Eaton <jwe@octave.org>
parents:
13983
diff
changeset
|
453 { |
99f039289e95
also clean up float fftw planner
John W. Eaton <jwe@octave.org>
parents:
13983
diff
changeset
|
454 instance = new octave_float_fftw_planner (); |
99f039289e95
also clean up float fftw planner
John W. Eaton <jwe@octave.org>
parents:
13983
diff
changeset
|
455 |
99f039289e95
also clean up float fftw planner
John W. Eaton <jwe@octave.org>
parents:
13983
diff
changeset
|
456 if (instance) |
99f039289e95
also clean up float fftw planner
John W. Eaton <jwe@octave.org>
parents:
13983
diff
changeset
|
457 singleton_cleanup_list::add (cleanup_instance); |
99f039289e95
also clean up float fftw planner
John W. Eaton <jwe@octave.org>
parents:
13983
diff
changeset
|
458 } |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
459 |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
460 if (! instance) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
461 { |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
462 (*current_liboctave_error_handler) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
463 ("unable to create octave_fftw_planner object!"); |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
464 |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
465 retval = false; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
466 } |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
467 |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
468 return retval; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
469 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
470 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
471 fftwf_plan |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
472 octave_float_fftw_planner::do_create_plan (int dir, const int rank, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
473 const dim_vector dims, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
474 octave_idx_type howmany, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
475 octave_idx_type stride, |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
476 octave_idx_type dist, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
477 const FloatComplex *in, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
478 FloatComplex *out) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
479 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
480 int which = (dir == FFTW_FORWARD) ? 0 : 1; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
481 fftwf_plan *cur_plan_p = &plan[which]; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
482 bool create_new_plan = false; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
483 bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
484 bool ioinplace = (in == out); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
485 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
486 // Don't create a new plan if we have a non SIMD plan already but |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
487 // can do SIMD. This prevents endlessly recreating plans if we |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
488 // change the alignment. |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
489 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
490 if (plan[which] == 0 || d[which] != dist || s[which] != stride |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
491 || r[which] != rank || h[which] != howmany |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
492 || ioinplace != inplace[which] |
20955
77f5591878bf
maint: Use '! expr' rather than '!expr' to conform to coding guidelines.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
493 || ((ioalign != simd_align[which]) ? ! ioalign : false)) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
494 create_new_plan = true; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
495 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
496 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
497 // We still might not have the same shape of array. |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
498 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
499 for (int i = 0; i < rank; i++) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
500 if (dims(i) != n[which](i)) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
501 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
502 create_new_plan = true; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
503 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
504 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
505 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
506 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
507 if (create_new_plan) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
508 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
509 d[which] = dist; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
510 s[which] = stride; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
511 r[which] = rank; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
512 h[which] = howmany; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
513 simd_align[which] = ioalign; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
514 inplace[which] = ioinplace; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
515 n[which] = dims; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
516 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
517 // Note reversal of dimensions for column major storage in FFTW. |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
518 octave_idx_type nn = 1; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
519 OCTAVE_LOCAL_BUFFER (int, tmp, rank); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
520 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
521 for (int i = 0, j = rank-1; i < rank; i++, j--) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
522 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
523 tmp[i] = dims(j); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
524 nn *= dims(j); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
525 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
526 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
527 int plan_flags = 0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
528 bool plan_destroys_in = true; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
529 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
530 switch (meth) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
531 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
532 case UNKNOWN: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
533 case ESTIMATE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
534 plan_flags |= FFTW_ESTIMATE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
535 plan_destroys_in = false; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
536 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
537 case MEASURE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
538 plan_flags |= FFTW_MEASURE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
539 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
540 case PATIENT: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
541 plan_flags |= FFTW_PATIENT; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
542 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
543 case EXHAUSTIVE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
544 plan_flags |= FFTW_EXHAUSTIVE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
545 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
546 case HYBRID: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
547 if (nn < 8193) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
548 plan_flags |= FFTW_MEASURE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
549 else |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
550 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
551 plan_flags |= FFTW_ESTIMATE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
552 plan_destroys_in = false; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
553 } |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
554 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
555 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
556 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
557 if (ioalign) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
558 plan_flags &= ~FFTW_UNALIGNED; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
559 else |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
560 plan_flags |= FFTW_UNALIGNED; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
561 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
562 if (*cur_plan_p) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
563 fftwf_destroy_plan (*cur_plan_p); |
3828 | 564 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
565 if (plan_destroys_in) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
566 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
567 // Create matrix with the same size and 16-byte alignment as input |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
568 OCTAVE_LOCAL_BUFFER (FloatComplex, itmp, nn * howmany + 32); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
569 itmp = reinterpret_cast<FloatComplex *> |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
570 (((reinterpret_cast<ptrdiff_t>(itmp) + 15) & ~ 0xF) + |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
571 ((reinterpret_cast<ptrdiff_t> (in)) & 0xF)); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
572 |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
573 *cur_plan_p = |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
574 fftwf_plan_many_dft (rank, tmp, howmany, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
575 reinterpret_cast<fftwf_complex *> (itmp), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
576 0, stride, dist, |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
577 reinterpret_cast<fftwf_complex *> (out), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
578 0, stride, dist, dir, plan_flags); |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
579 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
580 else |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
581 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
582 *cur_plan_p = |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
583 fftwf_plan_many_dft (rank, tmp, howmany, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
584 reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *> (in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
585 0, stride, dist, |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
586 reinterpret_cast<fftwf_complex *> (out), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
587 0, stride, dist, dir, plan_flags); |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
588 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
589 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
590 if (*cur_plan_p == 0) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
591 (*current_liboctave_error_handler) ("Error creating fftw plan"); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
592 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
593 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
594 return *cur_plan_p; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
595 } |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
596 |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
597 fftwf_plan |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
598 octave_float_fftw_planner::do_create_plan (const int rank, |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
599 const dim_vector dims, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
600 octave_idx_type howmany, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
601 octave_idx_type stride, |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
602 octave_idx_type dist, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
603 const float *in, FloatComplex *out) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
604 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
605 fftwf_plan *cur_plan_p = &rplan; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
606 bool create_new_plan = false; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
607 bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
608 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
609 // Don't create a new plan if we have a non SIMD plan already but |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
610 // can do SIMD. This prevents endlessly recreating plans if we |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
611 // change the alignment. |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
612 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
613 if (rplan == 0 || rd != dist || rs != stride || rr != rank |
20955
77f5591878bf
maint: Use '! expr' rather than '!expr' to conform to coding guidelines.
Rik <rik@octave.org>
parents:
19697
diff
changeset
|
614 || rh != howmany || ((ioalign != rsimd_align) ? ! ioalign : false)) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
615 create_new_plan = true; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
616 else |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
617 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
618 // We still might not have the same shape of array. |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
619 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
620 for (int i = 0; i < rank; i++) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
621 if (dims(i) != rn(i)) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
622 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
623 create_new_plan = true; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
624 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
625 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
626 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
627 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
628 if (create_new_plan) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
629 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
630 rd = dist; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
631 rs = stride; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
632 rr = rank; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
633 rh = howmany; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
634 rsimd_align = ioalign; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
635 rn = dims; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
636 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
637 // Note reversal of dimensions for column major storage in FFTW. |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
638 octave_idx_type nn = 1; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
639 OCTAVE_LOCAL_BUFFER (int, tmp, rank); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
640 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
641 for (int i = 0, j = rank-1; i < rank; i++, j--) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
642 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
643 tmp[i] = dims(j); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
644 nn *= dims(j); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
645 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
646 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
647 int plan_flags = 0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
648 bool plan_destroys_in = true; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
649 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
650 switch (meth) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
651 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
652 case UNKNOWN: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
653 case ESTIMATE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
654 plan_flags |= FFTW_ESTIMATE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
655 plan_destroys_in = false; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
656 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
657 case MEASURE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
658 plan_flags |= FFTW_MEASURE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
659 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
660 case PATIENT: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
661 plan_flags |= FFTW_PATIENT; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
662 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
663 case EXHAUSTIVE: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
664 plan_flags |= FFTW_EXHAUSTIVE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
665 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
666 case HYBRID: |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
667 if (nn < 8193) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
668 plan_flags |= FFTW_MEASURE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
669 else |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
670 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
671 plan_flags |= FFTW_ESTIMATE; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
672 plan_destroys_in = false; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
673 } |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
674 break; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
675 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
676 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
677 if (ioalign) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
678 plan_flags &= ~FFTW_UNALIGNED; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
679 else |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
680 plan_flags |= FFTW_UNALIGNED; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
681 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
682 if (*cur_plan_p) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
683 fftwf_destroy_plan (*cur_plan_p); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
684 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
685 if (plan_destroys_in) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
686 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
687 // Create matrix with the same size and 16-byte alignment as input |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
688 OCTAVE_LOCAL_BUFFER (float, itmp, nn + 32); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
689 itmp = reinterpret_cast<float *> |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
690 (((reinterpret_cast<ptrdiff_t>(itmp) + 15) & ~ 0xF) + |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
691 ((reinterpret_cast<ptrdiff_t> (in)) & 0xF)); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
692 |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
693 *cur_plan_p = |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
694 fftwf_plan_many_dft_r2c (rank, tmp, howmany, itmp, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
695 0, stride, dist, |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
696 reinterpret_cast<fftwf_complex *> (out), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
697 0, stride, dist, plan_flags); |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
698 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
699 else |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
700 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
701 *cur_plan_p = |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
702 fftwf_plan_many_dft_r2c (rank, tmp, howmany, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
703 (const_cast<float *> (in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
704 0, stride, dist, |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
705 reinterpret_cast<fftwf_complex *> (out), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
706 0, stride, dist, plan_flags); |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
707 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
708 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
709 if (*cur_plan_p == 0) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
710 (*current_liboctave_error_handler) ("Error creating fftw plan"); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
711 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
712 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
713 return *cur_plan_p; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
714 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
715 |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
716 octave_float_fftw_planner::FftwMethod |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
717 octave_float_fftw_planner::do_method (void) |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
718 { |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
719 return meth; |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
720 } |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
721 |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
722 octave_float_fftw_planner::FftwMethod |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
723 octave_float_fftw_planner::do_method (FftwMethod _meth) |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
724 { |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
725 FftwMethod ret = meth; |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
726 if (_meth == ESTIMATE || _meth == MEASURE |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
727 || _meth == PATIENT || _meth == EXHAUSTIVE |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
728 || _meth == HYBRID) |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
729 { |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
730 if (meth != _meth) |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
731 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
732 meth = _meth; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
733 if (rplan) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
734 fftwf_destroy_plan (rplan); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
735 if (plan[0]) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
736 fftwf_destroy_plan (plan[0]); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
737 if (plan[1]) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
738 fftwf_destroy_plan (plan[1]); |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
739 rplan = plan[0] = plan[1] = 0; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
740 } |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
741 } |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
742 else |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
743 ret = UNKNOWN; |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
744 return ret; |
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
745 } |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
746 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
747 template <class T> |
4775 | 748 static inline void |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
749 convert_packcomplex_1d (T *out, size_t nr, size_t nc, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
750 octave_idx_type stride, octave_idx_type dist) |
4773 | 751 { |
10142
829e69ec3110
make OCTAVE_QUIT a function
Jaroslav Hajek <highegg@gmail.com>
parents:
9523
diff
changeset
|
752 octave_quit (); |
4785 | 753 |
754 // Fill in the missing data. | |
755 | |
4773 | 756 for (size_t i = 0; i < nr; i++) |
757 for (size_t j = nc/2+1; j < nc; j++) | |
15018
3d8ace26c5b4
maint: Use Octave coding conventions for cuddled parentheses in liboctave/.
Rik <rik@octave.org>
parents:
14138
diff
changeset
|
758 out[j*stride + i*dist] = conj (out[(nc - j)*stride + i*dist]); |
4785 | 759 |
10142
829e69ec3110
make OCTAVE_QUIT a function
Jaroslav Hajek <highegg@gmail.com>
parents:
9523
diff
changeset
|
760 octave_quit (); |
4773 | 761 } |
762 | |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
763 template <class T> |
4775 | 764 static inline void |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
765 convert_packcomplex_Nd (T *out, const dim_vector &dv) |
3828 | 766 { |
4773 | 767 size_t nc = dv(0); |
768 size_t nr = dv(1); | |
4808 | 769 size_t np = (dv.length () > 2 ? dv.numel () / nc / nr : 1); |
4773 | 770 size_t nrp = nr * np; |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
771 T *ptr1, *ptr2; |
4773 | 772 |
10142
829e69ec3110
make OCTAVE_QUIT a function
Jaroslav Hajek <highegg@gmail.com>
parents:
9523
diff
changeset
|
773 octave_quit (); |
4785 | 774 |
775 // Create space for the missing elements. | |
776 | |
4773 | 777 for (size_t i = 0; i < nrp; i++) |
778 { | |
779 ptr1 = out + i * (nc/2 + 1) + nrp*((nc-1)/2); | |
780 ptr2 = out + i * nc; | |
781 for (size_t j = 0; j < nc/2+1; j++) | |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
782 *ptr2++ = *ptr1++; |
4773 | 783 } |
784 | |
10142
829e69ec3110
make OCTAVE_QUIT a function
Jaroslav Hajek <highegg@gmail.com>
parents:
9523
diff
changeset
|
785 octave_quit (); |
4785 | 786 |
787 // Fill in the missing data for the rank = 2 case directly for speed. | |
788 | |
4773 | 789 for (size_t i = 0; i < np; i++) |
790 { | |
791 for (size_t j = 1; j < nr; j++) | |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
792 for (size_t k = nc/2+1; k < nc; k++) |
15018
3d8ace26c5b4
maint: Use Octave coding conventions for cuddled parentheses in liboctave/.
Rik <rik@octave.org>
parents:
14138
diff
changeset
|
793 out[k + (j + i*nr)*nc] = conj (out[nc - k + ((i+1)*nr - j)*nc]); |
4773 | 794 |
795 for (size_t j = nc/2+1; j < nc; j++) | |
15018
3d8ace26c5b4
maint: Use Octave coding conventions for cuddled parentheses in liboctave/.
Rik <rik@octave.org>
parents:
14138
diff
changeset
|
796 out[j + i*nr*nc] = conj (out[(i*nr+1)*nc - j]); |
4773 | 797 } |
798 | |
10142
829e69ec3110
make OCTAVE_QUIT a function
Jaroslav Hajek <highegg@gmail.com>
parents:
9523
diff
changeset
|
799 octave_quit (); |
4785 | 800 |
801 // Now do the permutations needed for rank > 2 cases. | |
802 | |
4773 | 803 size_t jstart = dv(0) * dv(1); |
804 size_t kstep = dv(0); | |
805 size_t nel = dv.numel (); | |
4785 | 806 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
807 for (int inner = 2; inner < dv.length (); inner++) |
4773 | 808 { |
809 size_t jmax = jstart * dv(inner); | |
810 for (size_t i = 0; i < nel; i+=jmax) | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
811 for (size_t j = jstart, jj = jmax-jstart; j < jj; |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
812 j+=jstart, jj-=jstart) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
813 for (size_t k = 0; k < jstart; k+= kstep) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
814 for (size_t l = nc/2+1; l < nc; l++) |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
815 { |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
816 T tmp = out[i+ j + k + l]; |
15018
3d8ace26c5b4
maint: Use Octave coding conventions for cuddled parentheses in liboctave/.
Rik <rik@octave.org>
parents:
14138
diff
changeset
|
817 out[i + j + k + l] = out[i + jj + k + l]; |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
818 out[i + jj + k + l] = tmp; |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
819 } |
4773 | 820 jstart = jmax; |
821 } | |
4785 | 822 |
10142
829e69ec3110
make OCTAVE_QUIT a function
Jaroslav Hajek <highegg@gmail.com>
parents:
9523
diff
changeset
|
823 octave_quit (); |
4773 | 824 } |
825 | |
826 int | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
827 octave_fftw::fft (const double *in, Complex *out, size_t npts, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
828 size_t nsamples, octave_idx_type stride, octave_idx_type dist) |
4773 | 829 { |
830 dist = (dist < 0 ? npts : dist); | |
831 | |
10350
12884915a8e4
merge MArray classes & improve Array interface
Jaroslav Hajek <highegg@gmail.com>
parents:
10314
diff
changeset
|
832 dim_vector dv (npts, 1); |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
833 fftw_plan plan = octave_fftw_planner::create_plan (1, dv, nsamples, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
834 stride, dist, in, out); |
4773 | 835 |
836 fftw_execute_dft_r2c (plan, (const_cast<double *>(in)), | |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
837 reinterpret_cast<fftw_complex *> (out)); |
4773 | 838 |
4809 | 839 // Need to create other half of the transform. |
840 | |
4773 | 841 convert_packcomplex_1d (out, nsamples, npts, stride, dist); |
3828 | 842 |
843 return 0; | |
844 } | |
845 | |
846 int | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
847 octave_fftw::fft (const Complex *in, Complex *out, size_t npts, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
848 size_t nsamples, octave_idx_type stride, octave_idx_type dist) |
3828 | 849 { |
4773 | 850 dist = (dist < 0 ? npts : dist); |
851 | |
10350
12884915a8e4
merge MArray classes & improve Array interface
Jaroslav Hajek <highegg@gmail.com>
parents:
10314
diff
changeset
|
852 dim_vector dv (npts, 1); |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
853 fftw_plan plan = octave_fftw_planner::create_plan (FFTW_FORWARD, 1, dv, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
854 nsamples, stride, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
855 dist, in, out); |
4773 | 856 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
857 fftw_execute_dft (plan, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
858 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
859 reinterpret_cast<fftw_complex *> (out)); |
4773 | 860 |
861 return 0; | |
862 } | |
863 | |
864 int | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
865 octave_fftw::ifft (const Complex *in, Complex *out, size_t npts, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
866 size_t nsamples, octave_idx_type stride, |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
867 octave_idx_type dist) |
4773 | 868 { |
869 dist = (dist < 0 ? npts : dist); | |
870 | |
10350
12884915a8e4
merge MArray classes & improve Array interface
Jaroslav Hajek <highegg@gmail.com>
parents:
10314
diff
changeset
|
871 dim_vector dv (npts, 1); |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
872 fftw_plan plan = octave_fftw_planner::create_plan (FFTW_BACKWARD, 1, dv, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
873 nsamples, stride, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
874 dist, in, out); |
4773 | 875 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
876 fftw_execute_dft (plan, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
877 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
878 reinterpret_cast<fftw_complex *> (out)); |
3828 | 879 |
880 const Complex scale = npts; | |
4773 | 881 for (size_t j = 0; j < nsamples; j++) |
882 for (size_t i = 0; i < npts; i++) | |
883 out[i*stride + j*dist] /= scale; | |
3828 | 884 |
885 return 0; | |
886 } | |
887 | |
888 int | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
889 octave_fftw::fftNd (const double *in, Complex *out, const int rank, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
890 const dim_vector &dv) |
3828 | 891 { |
5275 | 892 octave_idx_type dist = 1; |
4773 | 893 for (int i = 0; i < rank; i++) |
894 dist *= dv(i); | |
895 | |
896 // Fool with the position of the start of the output matrix, so that | |
4809 | 897 // creating other half of the matrix won't cause cache problems. |
898 | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
899 octave_idx_type offset = (dv.numel () / dv(0)) * ((dv(0) - 1) / 2); |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
900 |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
901 fftw_plan plan = octave_fftw_planner::create_plan (rank, dv, 1, 1, dist, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
902 in, out + offset); |
4773 | 903 |
904 fftw_execute_dft_r2c (plan, (const_cast<double *>(in)), | |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
905 reinterpret_cast<fftw_complex *> (out+ offset)); |
4773 | 906 |
4809 | 907 // Need to create other half of the transform. |
908 | |
4773 | 909 convert_packcomplex_Nd (out, dv); |
3828 | 910 |
911 return 0; | |
912 } | |
913 | |
914 int | |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
915 octave_fftw::fftNd (const Complex *in, Complex *out, const int rank, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
916 const dim_vector &dv) |
3828 | 917 { |
5275 | 918 octave_idx_type dist = 1; |
4773 | 919 for (int i = 0; i < rank; i++) |
920 dist *= dv(i); | |
921 | |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
922 fftw_plan plan = octave_fftw_planner::create_plan (FFTW_FORWARD, rank, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
923 dv, 1, 1, dist, in, out); |
4773 | 924 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
925 fftw_execute_dft (plan, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
926 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
927 reinterpret_cast<fftw_complex *> (out)); |
4773 | 928 |
929 return 0; | |
930 } | |
3828 | 931 |
4773 | 932 int |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
933 octave_fftw::ifftNd (const Complex *in, Complex *out, const int rank, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
934 const dim_vector &dv) |
4773 | 935 { |
5275 | 936 octave_idx_type dist = 1; |
4773 | 937 for (int i = 0; i < rank; i++) |
938 dist *= dv(i); | |
939 | |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
940 fftw_plan plan = octave_fftw_planner::create_plan (FFTW_BACKWARD, rank, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
941 dv, 1, 1, dist, in, out); |
4773 | 942 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
943 fftw_execute_dft (plan, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
944 reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
945 reinterpret_cast<fftw_complex *> (out)); |
4773 | 946 |
947 const size_t npts = dv.numel (); | |
3828 | 948 const Complex scale = npts; |
949 for (size_t i = 0; i < npts; i++) | |
4773 | 950 out[i] /= scale; |
3828 | 951 |
952 return 0; | |
953 } | |
954 | |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
955 int |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
956 octave_fftw::fft (const float *in, FloatComplex *out, size_t npts, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
957 size_t nsamples, octave_idx_type stride, octave_idx_type dist) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
958 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
959 dist = (dist < 0 ? npts : dist); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
960 |
10350
12884915a8e4
merge MArray classes & improve Array interface
Jaroslav Hajek <highegg@gmail.com>
parents:
10314
diff
changeset
|
961 dim_vector dv (npts, 1); |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
962 fftwf_plan plan = octave_float_fftw_planner::create_plan (1, dv, nsamples, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
963 stride, dist, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
964 in, out); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
965 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
966 fftwf_execute_dft_r2c (plan, (const_cast<float *>(in)), |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
967 reinterpret_cast<fftwf_complex *> (out)); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
968 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
969 // Need to create other half of the transform. |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
970 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
971 convert_packcomplex_1d (out, nsamples, npts, stride, dist); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
972 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
973 return 0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
974 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
975 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
976 int |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
977 octave_fftw::fft (const FloatComplex *in, FloatComplex *out, size_t npts, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
978 size_t nsamples, octave_idx_type stride, octave_idx_type dist) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
979 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
980 dist = (dist < 0 ? npts : dist); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
981 |
10350
12884915a8e4
merge MArray classes & improve Array interface
Jaroslav Hajek <highegg@gmail.com>
parents:
10314
diff
changeset
|
982 dim_vector dv (npts, 1); |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
983 fftwf_plan plan = octave_float_fftw_planner::create_plan (FFTW_FORWARD, 1, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
984 dv, nsamples, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
985 stride, dist, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
986 in, out); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
987 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
988 fftwf_execute_dft (plan, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
989 reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
990 reinterpret_cast<fftwf_complex *> (out)); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
991 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
992 return 0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
993 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
994 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
995 int |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
996 octave_fftw::ifft (const FloatComplex *in, FloatComplex *out, size_t npts, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
997 size_t nsamples, octave_idx_type stride, |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
998 octave_idx_type dist) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
999 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1000 dist = (dist < 0 ? npts : dist); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1001 |
10350
12884915a8e4
merge MArray classes & improve Array interface
Jaroslav Hajek <highegg@gmail.com>
parents:
10314
diff
changeset
|
1002 dim_vector dv (npts, 1); |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
1003 fftwf_plan plan = octave_float_fftw_planner::create_plan (FFTW_BACKWARD, 1, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1004 dv, nsamples, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1005 stride, dist, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1006 in, out); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1007 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
1008 fftwf_execute_dft (plan, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1009 reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1010 reinterpret_cast<fftwf_complex *> (out)); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1011 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1012 const FloatComplex scale = npts; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1013 for (size_t j = 0; j < nsamples; j++) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1014 for (size_t i = 0; i < npts; i++) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1015 out[i*stride + j*dist] /= scale; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1016 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1017 return 0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1018 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1019 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1020 int |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
1021 octave_fftw::fftNd (const float *in, FloatComplex *out, const int rank, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1022 const dim_vector &dv) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1023 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1024 octave_idx_type dist = 1; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1025 for (int i = 0; i < rank; i++) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1026 dist *= dv(i); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1027 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1028 // Fool with the position of the start of the output matrix, so that |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1029 // creating other half of the matrix won't cause cache problems. |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1030 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
1031 octave_idx_type offset = (dv.numel () / dv(0)) * ((dv(0) - 1) / 2); |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
1032 |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
1033 fftwf_plan plan = octave_float_fftw_planner::create_plan (rank, dv, 1, 1, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1034 dist, in, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1035 out + offset); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1036 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1037 fftwf_execute_dft_r2c (plan, (const_cast<float *>(in)), |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1038 reinterpret_cast<fftwf_complex *> (out+ offset)); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1039 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1040 // Need to create other half of the transform. |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1041 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1042 convert_packcomplex_Nd (out, dv); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1043 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1044 return 0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1045 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1046 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1047 int |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
1048 octave_fftw::fftNd (const FloatComplex *in, FloatComplex *out, const int rank, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1049 const dim_vector &dv) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1050 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1051 octave_idx_type dist = 1; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1052 for (int i = 0; i < rank; i++) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1053 dist *= dv(i); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1054 |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
1055 fftwf_plan plan = octave_float_fftw_planner::create_plan (FFTW_FORWARD, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1056 rank, dv, 1, 1, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1057 dist, in, out); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1058 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
1059 fftwf_execute_dft (plan, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1060 reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1061 reinterpret_cast<fftwf_complex *> (out)); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1062 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1063 return 0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1064 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1065 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1066 int |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
1067 octave_fftw::ifftNd (const FloatComplex *in, FloatComplex *out, const int rank, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1068 const dim_vector &dv) |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1069 { |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1070 octave_idx_type dist = 1; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1071 for (int i = 0; i < rank; i++) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1072 dist *= dv(i); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1073 |
9516
fb933db0c517
convert fftw planner classes to singleton objects
John W. Eaton <jwe@octave.org>
parents:
8920
diff
changeset
|
1074 fftwf_plan plan = octave_float_fftw_planner::create_plan (FFTW_BACKWARD, |
10314
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1075 rank, dv, 1, 1, |
07ebe522dac2
untabify liboctave C++ sources
John W. Eaton <jwe@octave.org>
parents:
10158
diff
changeset
|
1076 dist, in, out); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1077 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
1078 fftwf_execute_dft (plan, |
17769
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1079 reinterpret_cast<fftwf_complex *> (const_cast<FloatComplex *>(in)), |
49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
1080 reinterpret_cast<fftwf_complex *> (out)); |
7789
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1081 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1082 const size_t npts = dv.numel (); |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1083 const FloatComplex scale = npts; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1084 for (size_t i = 0; i < npts; i++) |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1085 out[i] /= scale; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1086 |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1087 return 0; |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1088 } |
82be108cc558
First attempt at single precision tyeps
David Bateman <dbateman@free.fr>
parents:
7017
diff
changeset
|
1089 |
3828 | 1090 #endif |