1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486
|
/*
This file is part of Octave.
Octave is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2, or (at your option) any
later version.
Octave is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with Octave; see the file COPYING. If not, write to the Free
Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
02110-1301, USA.
*/
#ifdef HAVE_CONFIG_H
#include <config.h>
#endif
#if defined (HAVE_FFTW3)
#include <iostream>
#include <vector>
#include "lo-error.h"
#include "oct-fftw.h"
#include "quit.h"
// Helper class to create and cache fftw plans for both 1d and
// 2d. This implementation uses FFTW_ESTIMATE to create the plans,
// which in theory is suboptimal, but provides quite reasonable
// performance.
// Also note that if FFTW_ESTIMATE is not used the planner in FFTW3
// destroys the input and output arrays. So with the form of the
// current code we definitely want FFTW_ESTIMATE!! However, we use any
// wsidom that is available, either in a FFTW3 system wide file or as
// supplied by the user.
// XXX FIXME XXX -- if we can ensure 16 byte alignment in Array<T>
// (<T> *data) the FFTW3 can use SIMD instructions for further
// acceleration.
// Note that it is profitable to store the FFTW3 plans, for small
// ffts.
class
octave_fftw_planner
{
public:
octave_fftw_planner (void);
fftw_plan create_plan (int dir, const int rank, const dim_vector dims,
int howmany, int stride, int dist,
const Complex *in, Complex *out);
fftw_plan create_plan (const int rank, const dim_vector dims,
int howmany, int stride, int dist,
const double *in, Complex *out);
private:
int plan_flags;
// XXX FIXME XXX -- perhaps this should be split into two classes?
// Plan for fft and ifft of complex values
fftw_plan plan[2];
// dist
int d[2];
// stride
int s[2];
// rank
int r[2];
// howmany
int h[2];
// dims
dim_vector n[2];
bool simd_align[2];
bool inplace[2];
// Plan for fft of real values
fftw_plan rplan;
// dist
int rd;
// stride
int rs;
// rank
int rr;
// howmany
int rh;
// dims
dim_vector rn;
bool rsimd_align;
};
octave_fftw_planner::octave_fftw_planner (void)
{
plan_flags = FFTW_ESTIMATE;
plan[0] = plan[1] = 0;
d[0] = d[1] = s[0] = s[1] = r[0] = r[1] = h[0] = h[1] = 0;
simd_align[0] = simd_align[1] = false;
inplace[0] = inplace[1] = false;
n[0] = n[1] = dim_vector ();
rplan = 0;
rd = rs = rr = rh = 0;
rsimd_align = false;
rn = dim_vector ();
// If we have a system wide wisdom file, import it.
fftw_import_system_wisdom ();
}
#define CHECK_SIMD_ALIGNMENT(x) \
((reinterpret_cast<ptrdiff_t> (x)) & 0xF == 0)
fftw_plan
octave_fftw_planner::create_plan (int dir, const int rank,
const dim_vector dims, int howmany,
int stride, int dist,
const Complex *in, Complex *out)
{
int which = (dir == FFTW_FORWARD) ? 0 : 1;
fftw_plan *cur_plan_p = &plan[which];
bool create_new_plan = false;
bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out);
bool ioinplace = (in == out);
// Don't create a new plan if we have a non SIMD plan already but
// can do SIMD. This prevents endlessly recreating plans if we
// change the alignment.
if (plan[which] == 0 || d[which] != dist || s[which] != stride
|| r[which] != rank || h[which] != howmany
|| ioinplace != inplace[which]
|| ((ioalign != simd_align[which]) ? !ioalign : false))
create_new_plan = true;
else
{
// We still might not have the same shape of array.
for (int i = 0; i < rank; i++)
if (dims(i) != n[which](i))
{
create_new_plan = true;
break;
}
}
if (create_new_plan)
{
d[which] = dist;
s[which] = stride;
r[which] = rank;
h[which] = howmany;
simd_align[which] = ioalign;
inplace[which] = ioinplace;
n[which] = dims;
if (ioalign)
plan_flags &= ~FFTW_UNALIGNED;
else
plan_flags |= FFTW_UNALIGNED;
if (*cur_plan_p)
fftw_destroy_plan (*cur_plan_p);
// Note reversal of dimensions for column major storage in FFTW.
OCTAVE_LOCAL_BUFFER (int, tmp, rank);
for (int i = 0, j = rank-1; i < rank; i++, j--)
tmp[i] = dims(j);
*cur_plan_p =
fftw_plan_many_dft (rank, tmp, howmany,
reinterpret_cast<fftw_complex *> (const_cast<Complex *> (in)),
0, stride, dist, reinterpret_cast<fftw_complex *> (out),
0, stride, dist, dir, plan_flags);
if (*cur_plan_p == 0)
(*current_liboctave_error_handler) ("Error creating fftw plan");
}
return *cur_plan_p;
}
fftw_plan
octave_fftw_planner::create_plan (const int rank, const dim_vector dims,
int howmany, int stride, int dist,
const double *in, Complex *out)
{
fftw_plan *cur_plan_p = &rplan;
bool create_new_plan = false;
bool ioalign = CHECK_SIMD_ALIGNMENT (in) && CHECK_SIMD_ALIGNMENT (out);
// Don't create a new plan if we have a non SIMD plan already but
// can do SIMD. This prevents endlessly recreating plans if we
// change the alignment.
if (rplan == 0 || rd != dist || rs != stride || rr != rank
|| rh != howmany || ((ioalign != rsimd_align) ? !ioalign : false))
create_new_plan = true;
else
{
// We still might not have the same shape of array.
for (int i = 0; i < rank; i++)
if (dims(i) != rn(i))
{
create_new_plan = true;
break;
}
}
if (create_new_plan)
{
rd = dist;
rs = stride;
rr = rank;
rh = howmany;
rsimd_align = ioalign;
rn = dims;
if (ioalign)
plan_flags &= ~FFTW_UNALIGNED;
else
plan_flags |= FFTW_UNALIGNED;
if (*cur_plan_p)
fftw_destroy_plan (*cur_plan_p);
// Note reversal of dimensions for column major storage in FFTW.
OCTAVE_LOCAL_BUFFER (int, tmp, rank);
for (int i = 0, j = rank-1; i < rank; i++, j--)
tmp[i] = dims(j);
*cur_plan_p =
fftw_plan_many_dft_r2c (rank, tmp, howmany,
(const_cast<double *> (in)),
0, stride, dist, reinterpret_cast<fftw_complex *> (out),
0, stride, dist, plan_flags);
if (*cur_plan_p == 0)
(*current_liboctave_error_handler) ("Error creating fftw plan");
}
return *cur_plan_p;
}
static octave_fftw_planner fftw_planner;
static inline void
convert_packcomplex_1d (Complex *out, size_t nr, size_t nc,
int stride, int dist)
{
OCTAVE_QUIT;
// Fill in the missing data.
for (size_t i = 0; i < nr; i++)
for (size_t j = nc/2+1; j < nc; j++)
out[j*stride + i*dist] = conj(out[(nc - j)*stride + i*dist]);
OCTAVE_QUIT;
}
static inline void
convert_packcomplex_Nd (Complex *out, const dim_vector &dv)
{
size_t nc = dv(0);
size_t nr = dv(1);
size_t np = (dv.length () > 2 ? dv.numel () / nc / nr : 1);
size_t nrp = nr * np;
Complex *ptr1, *ptr2;
OCTAVE_QUIT;
// Create space for the missing elements.
for (size_t i = 0; i < nrp; i++)
{
ptr1 = out + i * (nc/2 + 1) + nrp*((nc-1)/2);
ptr2 = out + i * nc;
for (size_t j = 0; j < nc/2+1; j++)
*ptr2++ = *ptr1++;
}
OCTAVE_QUIT;
// Fill in the missing data for the rank = 2 case directly for speed.
for (size_t i = 0; i < np; i++)
{
for (size_t j = 1; j < nr; j++)
for (size_t k = nc/2+1; k < nc; k++)
out[k + (j + i*nr)*nc] = conj(out[nc - k + ((i+1)*nr - j)*nc]);
for (size_t j = nc/2+1; j < nc; j++)
out[j + i*nr*nc] = conj(out[(i*nr+1)*nc - j]);
}
OCTAVE_QUIT;
// Now do the permutations needed for rank > 2 cases.
size_t jstart = dv(0) * dv(1);
size_t kstep = dv(0);
size_t nel = dv.numel ();
for (int inner = 2; inner < dv.length (); inner++)
{
size_t jmax = jstart * dv(inner);
for (size_t i = 0; i < nel; i+=jmax)
for (size_t j = jstart, jj = jmax-jstart; j < jj;
j+=jstart, jj-=jstart)
for (size_t k = 0; k < jstart; k+= kstep)
for (size_t l = nc/2+1; l < nc; l++)
{
Complex tmp = out[i+ j + k + l];
out[i + j + k + l] = out[i + jj + k + l];
out[i + jj + k + l] = tmp;
}
jstart = jmax;
}
OCTAVE_QUIT;
}
int
octave_fftw::fft (const double *in, Complex *out, size_t npts,
size_t nsamples, int stride, int dist)
{
dist = (dist < 0 ? npts : dist);
dim_vector dv (npts);
fftw_plan plan = fftw_planner.create_plan (1, dv, nsamples, stride, dist,
in, out);
fftw_execute_dft_r2c (plan, (const_cast<double *>(in)),
reinterpret_cast<fftw_complex *> (out));
// Need to create other half of the transform.
convert_packcomplex_1d (out, nsamples, npts, stride, dist);
return 0;
}
int
octave_fftw::fft (const Complex *in, Complex *out, size_t npts,
size_t nsamples, int stride, int dist)
{
dist = (dist < 0 ? npts : dist);
dim_vector dv (npts);
fftw_plan plan = fftw_planner.create_plan (FFTW_FORWARD, 1, dv, nsamples,
stride, dist, in, out);
fftw_execute_dft (plan,
reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
reinterpret_cast<fftw_complex *> (out));
return 0;
}
int
octave_fftw::ifft (const Complex *in, Complex *out, size_t npts,
size_t nsamples, int stride, int dist)
{
dist = (dist < 0 ? npts : dist);
dim_vector dv (npts);
fftw_plan plan = fftw_planner.create_plan (FFTW_BACKWARD, 1, dv, nsamples,
stride, dist, in, out);
fftw_execute_dft (plan,
reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
reinterpret_cast<fftw_complex *> (out));
const Complex scale = npts;
for (size_t j = 0; j < nsamples; j++)
for (size_t i = 0; i < npts; i++)
out[i*stride + j*dist] /= scale;
return 0;
}
int
octave_fftw::fftNd (const double *in, Complex *out, const int rank,
const dim_vector &dv)
{
int dist = 1;
for (int i = 0; i < rank; i++)
dist *= dv(i);
// Fool with the position of the start of the output matrix, so that
// creating other half of the matrix won't cause cache problems.
int offset = (dv.numel () / dv(0)) * ((dv(0) - 1) / 2);
fftw_plan plan = fftw_planner.create_plan (rank, dv, 1, 1, dist,
in, out + offset);
fftw_execute_dft_r2c (plan, (const_cast<double *>(in)),
reinterpret_cast<fftw_complex *> (out+ offset));
// Need to create other half of the transform.
convert_packcomplex_Nd (out, dv);
return 0;
}
int
octave_fftw::fftNd (const Complex *in, Complex *out, const int rank,
const dim_vector &dv)
{
int dist = 1;
for (int i = 0; i < rank; i++)
dist *= dv(i);
fftw_plan plan = fftw_planner.create_plan (FFTW_FORWARD, rank, dv, 1, 1,
dist, in, out);
fftw_execute_dft (plan,
reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
reinterpret_cast<fftw_complex *> (out));
return 0;
}
int
octave_fftw::ifftNd (const Complex *in, Complex *out, const int rank,
const dim_vector &dv)
{
int dist = 1;
for (int i = 0; i < rank; i++)
dist *= dv(i);
fftw_plan plan = fftw_planner.create_plan (FFTW_BACKWARD, rank, dv, 1, 1,
dist, in, out);
fftw_execute_dft (plan,
reinterpret_cast<fftw_complex *> (const_cast<Complex *>(in)),
reinterpret_cast<fftw_complex *> (out));
const size_t npts = dv.numel ();
const Complex scale = npts;
for (size_t i = 0; i < npts; i++)
out[i] /= scale;
return 0;
}
#endif
/*
;;; Local Variables: ***
;;; mode: C++ ***
;;; End: ***
*/
|