source: src/spectral/fft.c @ 81abf91

feature/autosinkfeature/cnnfeature/cnn_orgfeature/constantqfeature/crepefeature/crepe_orgfeature/pitchshiftfeature/pydocstringsfeature/timestretchfix/ffmpeg5
Last change on this file since 81abf91 was b701179, checked in by Paul Brossier <piem@piem.org>, 7 years ago

src/*.c, wscript: remove trailing spaces

  • Property mode set to 100644
File size: 17.6 KB
RevLine 
[96fb8ad]1/*
[e6a78ea]2  Copyright (C) 2003-2009 Paul Brossier <piem@aubio.org>
[96fb8ad]3
[e6a78ea]4  This file is part of aubio.
[96fb8ad]5
[e6a78ea]6  aubio is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
[96fb8ad]10
[e6a78ea]11  aubio is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
14  GNU General Public License for more details.
15
16  You should have received a copy of the GNU General Public License
17  along with aubio.  If not, see <http://www.gnu.org/licenses/>.
[96fb8ad]18
19*/
20
21#include "aubio_priv.h"
[aadd27a]22#include "fvec.h"
23#include "cvec.h"
[96fb8ad]24#include "mathutils.h"
[32d6958]25#include "spectral/fft.h"
[96fb8ad]26
[54dd945]27#ifdef HAVE_FFTW3             // using FFTW3
[f3bee79]28/* note that <complex.h> is not included here but only in aubio_priv.h, so that
29 * c++ projects can still use their own complex definition. */
30#include <fftw3.h>
[d453a4a]31#include <pthread.h>
[f3bee79]32
33#ifdef HAVE_COMPLEX_H
[729a3c0]34#ifdef HAVE_FFTW3F
[f3bee79]35/** fft data type with complex.h and fftw3f */
36#define FFTW_TYPE fftwf_complex
37#else
38/** fft data type with complex.h and fftw3 */
39#define FFTW_TYPE fftw_complex
40#endif
41#else
[729a3c0]42#ifdef HAVE_FFTW3F
[f3bee79]43/** fft data type without complex.h and with fftw3f */
44#define FFTW_TYPE float
45#else
46/** fft data type without complex.h and with fftw */
47#define FFTW_TYPE double
48#endif
49#endif
50
51/** fft data type */
52typedef FFTW_TYPE fft_data_t;
53
[729a3c0]54#ifdef HAVE_FFTW3F
[a47cd35]55#define fftw_malloc            fftwf_malloc
56#define fftw_free              fftwf_free
57#define fftw_execute           fftwf_execute
58#define fftw_plan_dft_r2c_1d   fftwf_plan_dft_r2c_1d
59#define fftw_plan_dft_c2r_1d   fftwf_plan_dft_c2r_1d
60#define fftw_plan_r2r_1d       fftwf_plan_r2r_1d
61#define fftw_plan              fftwf_plan
62#define fftw_destroy_plan      fftwf_destroy_plan
[96fb8ad]63#endif
64
[729a3c0]65#ifdef HAVE_FFTW3F
[c204928]66#if HAVE_AUBIO_DOUBLE
[eeb7276]67#error "Using aubio in double precision with fftw3 in single precision"
[fd6b90f]68#endif /* HAVE_AUBIO_DOUBLE */
[152ed05]69#define real_t float
[eeb7276]70#elif defined (HAVE_FFTW3) /* HAVE_FFTW3F */
[c204928]71#if !HAVE_AUBIO_DOUBLE
[eeb7276]72#error "Using aubio in single precision with fftw3 in double precision"
[fd6b90f]73#endif /* HAVE_AUBIO_DOUBLE */
[152ed05]74#define real_t double
[fd6b90f]75#endif /* HAVE_FFTW3F */
[96fb8ad]76
[d453a4a]77// a global mutex for FFTW thread safety
78pthread_mutex_t aubio_fftw_mutex = PTHREAD_MUTEX_INITIALIZER;
79
[986131d]80#elif defined HAVE_ACCELERATE        // using ACCELERATE
[54dd945]81// https://developer.apple.com/library/mac/#documentation/Accelerate/Reference/vDSPRef/Reference/reference.html
82#include <Accelerate/Accelerate.h>
83
[39c4721]84#if !HAVE_AUBIO_DOUBLE
85#define aubio_vDSP_ctoz                vDSP_ctoz
86#define aubio_vDSP_fft_zrip            vDSP_fft_zrip
87#define aubio_vDSP_ztoc                vDSP_ztoc
88#define aubio_vDSP_zvmags              vDSP_zvmags
89#define aubio_vDSP_zvphas              vDSP_zvphas
90#define aubio_vDSP_vsadd               vDSP_vsadd
91#define aubio_vDSP_vsmul               vDSP_vsmul
92#define aubio_vDSP_create_fftsetup     vDSP_create_fftsetup
93#define aubio_vDSP_destroy_fftsetup    vDSP_destroy_fftsetup
94#define aubio_DSPComplex               DSPComplex
95#define aubio_DSPSplitComplex          DSPSplitComplex
96#define aubio_FFTSetup                 FFTSetup
97#define aubio_vvsqrt                   vvsqrtf
98#else
99#define aubio_vDSP_ctoz                vDSP_ctozD
100#define aubio_vDSP_fft_zrip            vDSP_fft_zripD
101#define aubio_vDSP_ztoc                vDSP_ztocD
102#define aubio_vDSP_zvmags              vDSP_zvmagsD
103#define aubio_vDSP_zvphas              vDSP_zvphasD
104#define aubio_vDSP_vsadd               vDSP_vsaddD
105#define aubio_vDSP_vsmul               vDSP_vsmulD
106#define aubio_vDSP_create_fftsetup     vDSP_create_fftsetupD
107#define aubio_vDSP_destroy_fftsetup    vDSP_destroy_fftsetupD
108#define aubio_DSPComplex               DSPDoubleComplex
109#define aubio_DSPSplitComplex          DSPDoubleSplitComplex
110#define aubio_FFTSetup                 FFTSetupD
111#define aubio_vvsqrt                   vvsqrt
112#endif /* HAVE_AUBIO_DOUBLE */
113
[986131d]114#elif defined HAVE_INTEL_IPP // using INTEL IPP
115
[0ad2e17]116#if !HAVE_AUBIO_DOUBLE
117#define aubio_IppFloat                 Ipp32f
118#define aubio_IppComplex               Ipp32fc
119#define aubio_FFTSpec                  FFTSpec_R_32f
120#define aubio_ippsMalloc_complex       ippsMalloc_32fc
121#define aubio_ippsFFTInit_R            ippsFFTInit_R_32f
122#define aubio_ippsFFTGetSize_R         ippsFFTGetSize_R_32f
123#define aubio_ippsFFTInv_CCSToR        ippsFFTInv_CCSToR_32f
124#define aubio_ippsFFTFwd_RToCCS        ippsFFTFwd_RToCCS_32f
[14e300e]125#define aubio_ippsAtan2                ippsAtan2_32f_A21
[0ad2e17]126#else /* HAVE_AUBIO_DOUBLE */
127#define aubio_IppFloat                 Ipp64f
128#define aubio_IppComplex               Ipp64fc
129#define aubio_FFTSpec                  FFTSpec_R_64f
130#define aubio_ippsMalloc_complex       ippsMalloc_64fc
131#define aubio_ippsFFTInit_R            ippsFFTInit_R_64f
132#define aubio_ippsFFTGetSize_R         ippsFFTGetSize_R_64f
133#define aubio_ippsFFTInv_CCSToR        ippsFFTInv_CCSToR_64f
134#define aubio_ippsFFTFwd_RToCCS        ippsFFTFwd_RToCCS_64f
[14e300e]135#define aubio_ippsAtan2                ippsAtan2_64f_A50
[0ad2e17]136#endif
137
[986131d]138
139#else // using OOURA
[54dd945]140// let's use ooura instead
[f50c9503]141extern void aubio_ooura_rdft(int, int, smpl_t *, int *, smpl_t *);
[54dd945]142
[986131d]143#endif
[729a3c0]144
[96fb8ad]145struct _aubio_fft_t {
[aadd27a]146  uint_t winsize;
147  uint_t fft_size;
[986131d]148
[54dd945]149#ifdef HAVE_FFTW3             // using FFTW3
[aadd27a]150  real_t *in, *out;
[4b6937b]151  fftw_plan pfw, pbw;
[986131d]152  fft_data_t * specdata; /* complex spectral data */
153
154#elif defined HAVE_ACCELERATE  // using ACCELERATE
[54dd945]155  int log2fftsize;
[39c4721]156  aubio_FFTSetup fftSetup;
157  aubio_DSPSplitComplex spec;
158  smpl_t *in, *out;
[b701179]159
[986131d]160#elif defined HAVE_INTEL_IPP  // using Intel IPP
161  smpl_t *in, *out;
162  Ipp8u* memSpec;
163  Ipp8u* memInit;
164  Ipp8u* memBuffer;
[0ad2e17]165  struct aubio_FFTSpec* fftSpec;
166  aubio_IppComplex* complexOut;
[54dd945]167#else                         // using OOURA
[5c6b264]168  smpl_t *in, *out;
169  smpl_t *w;
[729a3c0]170  int *ip;
[986131d]171#endif /* using OOURA */
172
[aadd27a]173  fvec_t * compspec;
[96fb8ad]174};
175
[729a3c0]176aubio_fft_t * new_aubio_fft (uint_t winsize) {
[a47cd35]177  aubio_fft_t * s = AUBIO_NEW(aubio_fft_t);
[d897aae]178  if ((sint_t)winsize < 2) {
179    AUBIO_ERR("fft: got winsize %d, but can not be < 2\n", winsize);
[a984461]180    goto beach;
181  }
[986131d]182
[729a3c0]183#ifdef HAVE_FFTW3
[8b3a7e7]184  uint_t i;
[aadd27a]185  s->winsize  = winsize;
[a47cd35]186  /* allocate memory */
[aadd27a]187  s->in       = AUBIO_ARRAY(real_t,winsize);
188  s->out      = AUBIO_ARRAY(real_t,winsize);
[d95ff38]189  s->compspec = new_fvec(winsize);
[a47cd35]190  /* create plans */
[d453a4a]191  pthread_mutex_lock(&aubio_fftw_mutex);
[237f632]192#ifdef HAVE_COMPLEX_H
[4b6937b]193  s->fft_size = winsize/2 + 1;
[a47cd35]194  s->specdata = (fft_data_t*)fftw_malloc(sizeof(fft_data_t)*s->fft_size);
[aadd27a]195  s->pfw = fftw_plan_dft_r2c_1d(winsize, s->in,  s->specdata, FFTW_ESTIMATE);
196  s->pbw = fftw_plan_dft_c2r_1d(winsize, s->specdata, s->out, FFTW_ESTIMATE);
[237f632]197#else
[aadd27a]198  s->fft_size = winsize;
[a47cd35]199  s->specdata = (fft_data_t*)fftw_malloc(sizeof(fft_data_t)*s->fft_size);
[aadd27a]200  s->pfw = fftw_plan_r2r_1d(winsize, s->in,  s->specdata, FFTW_R2HC, FFTW_ESTIMATE);
201  s->pbw = fftw_plan_r2r_1d(winsize, s->specdata, s->out, FFTW_HC2R, FFTW_ESTIMATE);
[237f632]202#endif
[d453a4a]203  pthread_mutex_unlock(&aubio_fftw_mutex);
[4f4299d]204  for (i = 0; i < s->winsize; i++) {
205    s->in[i] = 0.;
206    s->out[i] = 0.;
207  }
208  for (i = 0; i < s->fft_size; i++) {
209    s->specdata[i] = 0.;
210  }
[986131d]211
212#elif defined HAVE_ACCELERATE  // using ACCELERATE
[54dd945]213  s->winsize = winsize;
214  s->fft_size = winsize;
215  s->compspec = new_fvec(winsize);
[986131d]216  s->log2fftsize = aubio_power_of_two_order(s->fft_size);
[39c4721]217  s->in = AUBIO_ARRAY(smpl_t, s->fft_size);
218  s->out = AUBIO_ARRAY(smpl_t, s->fft_size);
219  s->spec.realp = AUBIO_ARRAY(smpl_t, s->fft_size/2);
220  s->spec.imagp = AUBIO_ARRAY(smpl_t, s->fft_size/2);
221  s->fftSetup = aubio_vDSP_create_fftsetup(s->log2fftsize, FFT_RADIX2);
[986131d]222
223#elif defined HAVE_INTEL_IPP  // using Intel IPP
224  const IppHintAlgorithm qualityHint = ippAlgHintAccurate; // OR ippAlgHintFast;
[b701179]225  const int flags = IPP_FFT_NODIV_BY_ANY; // we're scaling manually afterwards
[986131d]226  int order = aubio_power_of_two_order(winsize);
227  int sizeSpec, sizeInit, sizeBuffer;
228  IppStatus status;
229
230  if (winsize <= 4 || aubio_is_power_of_two(winsize) != 1)
231  {
232    AUBIO_ERR("intel IPP fft: can only create with sizes > 4 and power of two, requested %d,"
233      " try recompiling aubio with --enable-fftw3\n", winsize);
234    goto beach;
235  }
236
[0ad2e17]237  status = aubio_ippsFFTGetSize_R(order, flags, qualityHint,
[986131d]238      &sizeSpec, &sizeInit, &sizeBuffer);
239  if (status != ippStsNoErr) {
240    AUBIO_ERR("fft: failed to initialize fft. IPP error: %d\n", status);
241    goto beach;
242  }
243  s->fft_size = s->winsize = winsize;
244  s->compspec = new_fvec(winsize);
245  s->in = AUBIO_ARRAY(smpl_t, s->winsize);
246  s->out = AUBIO_ARRAY(smpl_t, s->winsize);
247  s->memSpec = ippsMalloc_8u(sizeSpec);
248  s->memBuffer = ippsMalloc_8u(sizeBuffer);
249  if (sizeInit > 0 ) {
250    s->memInit = ippsMalloc_8u(sizeInit);
251  }
[0ad2e17]252  s->complexOut = aubio_ippsMalloc_complex(s->fft_size / 2 + 1);
253  status = aubio_ippsFFTInit_R(
[986131d]254    &s->fftSpec, order, flags, qualityHint, s->memSpec, s->memInit);
255  if (status != ippStsNoErr) {
256    AUBIO_ERR("fft: failed to initialize. IPP error: %d\n", status);
257    goto beach;
258  }
259
[54dd945]260#else                         // using OOURA
[3c6f584]261  if (aubio_is_power_of_two(winsize) != 1) {
[1b57274]262    AUBIO_ERR("fft: can only create with sizes power of two, requested %d,"
263        " try recompiling aubio with --enable-fftw3\n", winsize);
[a984461]264    goto beach;
[3c6f584]265  }
[729a3c0]266  s->winsize = winsize;
267  s->fft_size = winsize / 2 + 1;
268  s->compspec = new_fvec(winsize);
[5c6b264]269  s->in    = AUBIO_ARRAY(smpl_t, s->winsize);
270  s->out   = AUBIO_ARRAY(smpl_t, s->winsize);
[729a3c0]271  s->ip    = AUBIO_ARRAY(int   , s->fft_size);
[5c6b264]272  s->w     = AUBIO_ARRAY(smpl_t, s->fft_size);
[729a3c0]273  s->ip[0] = 0;
[986131d]274#endif /* using OOURA */
275
[a47cd35]276  return s;
[986131d]277
[a984461]278beach:
279  AUBIO_FREE(s);
280  return NULL;
[96fb8ad]281}
282
283void del_aubio_fft(aubio_fft_t * s) {
[a47cd35]284  /* destroy data */
[54dd945]285#ifdef HAVE_FFTW3             // using FFTW3
[4b251ae]286  pthread_mutex_lock(&aubio_fftw_mutex);
[a47cd35]287  fftw_destroy_plan(s->pfw);
288  fftw_destroy_plan(s->pbw);
289  fftw_free(s->specdata);
[4b251ae]290  pthread_mutex_unlock(&aubio_fftw_mutex);
[986131d]291
292#elif defined HAVE_ACCELERATE // using ACCELERATE
[54dd945]293  AUBIO_FREE(s->spec.realp);
294  AUBIO_FREE(s->spec.imagp);
[39c4721]295  aubio_vDSP_destroy_fftsetup(s->fftSetup);
[986131d]296
297#elif defined HAVE_INTEL_IPP  // using Intel IPP
298  ippFree(s->memSpec);
299  ippFree(s->memInit);
300  ippFree(s->memBuffer);
301  ippFree(s->complexOut);
302
[54dd945]303#else                         // using OOURA
[729a3c0]304  AUBIO_FREE(s->w);
305  AUBIO_FREE(s->ip);
[986131d]306#endif
307
308  del_fvec(s->compspec);
[729a3c0]309  AUBIO_FREE(s->in);
[986131d]310  AUBIO_FREE(s->out);
[a47cd35]311  AUBIO_FREE(s);
[96fb8ad]312}
313
[feb694b]314void aubio_fft_do(aubio_fft_t * s, const fvec_t * input, cvec_t * spectrum) {
[aadd27a]315  aubio_fft_do_complex(s, input, s->compspec);
[14e300e]316  aubio_fft_get_spectrum(s->compspec, spectrum);
[96fb8ad]317}
318
[feb694b]319void aubio_fft_rdo(aubio_fft_t * s, const cvec_t * spectrum, fvec_t * output) {
[14e300e]320  aubio_fft_get_realimag(spectrum, s->compspec);
[aadd27a]321  aubio_fft_rdo_complex(s, s->compspec, output);
[96fb8ad]322}
323
[feb694b]324void aubio_fft_do_complex(aubio_fft_t * s, const fvec_t * input, fvec_t * compspec) {
[729a3c0]325  uint_t i;
[b5d32cb]326#ifndef HAVE_MEMCPY_HACKS
[729a3c0]327  for (i=0; i < s->winsize; i++) {
328    s->in[i] = input->data[i];
[d95ff38]329  }
[b5d32cb]330#else
331  memcpy(s->in, input->data, s->winsize * sizeof(smpl_t));
332#endif /* HAVE_MEMCPY_HACKS */
[986131d]333
[54dd945]334#ifdef HAVE_FFTW3             // using FFTW3
[d95ff38]335  fftw_execute(s->pfw);
[4b6937b]336#ifdef HAVE_COMPLEX_H
[d95ff38]337  compspec->data[0] = REAL(s->specdata[0]);
[729a3c0]338  for (i = 1; i < s->fft_size -1 ; i++) {
339    compspec->data[i] = REAL(s->specdata[i]);
340    compspec->data[compspec->length - i] = IMAG(s->specdata[i]);
[d95ff38]341  }
342  compspec->data[s->fft_size-1] = REAL(s->specdata[s->fft_size-1]);
[729a3c0]343#else /* HAVE_COMPLEX_H  */
344  for (i = 0; i < s->fft_size; i++) {
345    compspec->data[i] = s->specdata[i];
[237f632]346  }
[729a3c0]347#endif /* HAVE_COMPLEX_H */
[986131d]348
349#elif defined HAVE_ACCELERATE // using ACCELERATE
[54dd945]350  // convert real data to even/odd format used in vDSP
[39c4721]351  aubio_vDSP_ctoz((aubio_DSPComplex*)s->in, 2, &s->spec, 1, s->fft_size/2);
[54dd945]352  // compute the FFT
[39c4721]353  aubio_vDSP_fft_zrip(s->fftSetup, &s->spec, 1, s->log2fftsize, FFT_FORWARD);
[54dd945]354  // convert from vDSP complex split to [ r0, r1, ..., rN, iN-1, .., i2, i1]
355  compspec->data[0] = s->spec.realp[0];
356  compspec->data[s->fft_size / 2] = s->spec.imagp[0];
357  for (i = 1; i < s->fft_size / 2; i++) {
358    compspec->data[i] = s->spec.realp[i];
359    compspec->data[s->fft_size - i] = s->spec.imagp[i];
360  }
361  // apply scaling
362  smpl_t scale = 1./2.;
[39c4721]363  aubio_vDSP_vsmul(compspec->data, 1, &scale, compspec->data, 1, s->fft_size);
[986131d]364
365#elif defined HAVE_INTEL_IPP  // using Intel IPP
366
367  // apply fft
[0ad2e17]368  aubio_ippsFFTFwd_RToCCS(s->in, (aubio_IppFloat*)s->complexOut, s->fftSpec, s->memBuffer);
[986131d]369  // convert complex buffer to [ r0, r1, ..., rN, iN-1, .., i2, i1]
370  compspec->data[0] = s->complexOut[0].re;
371  compspec->data[s->fft_size / 2] = s->complexOut[s->fft_size / 2].re;
372  for (i = 1; i < s->fft_size / 2; i++) {
373    compspec->data[i] = s->complexOut[i].re;
374    compspec->data[s->fft_size - i] = s->complexOut[i].im;
375  }
376
[54dd945]377#else                         // using OOURA
[f50c9503]378  aubio_ooura_rdft(s->winsize, 1, s->in, s->ip, s->w);
[729a3c0]379  compspec->data[0] = s->in[0];
380  compspec->data[s->winsize / 2] = s->in[1];
381  for (i = 1; i < s->fft_size - 1; i++) {
382    compspec->data[i] = s->in[2 * i];
383    compspec->data[s->winsize - i] = - s->in[2 * i + 1];
384  }
[986131d]385#endif /* using OOURA */
[96fb8ad]386}
387
[feb694b]388void aubio_fft_rdo_complex(aubio_fft_t * s, const fvec_t * compspec, fvec_t * output) {
[729a3c0]389  uint_t i;
390#ifdef HAVE_FFTW3
[aadd27a]391  const smpl_t renorm = 1./(smpl_t)s->winsize;
[4b6937b]392#ifdef HAVE_COMPLEX_H
[d95ff38]393  s->specdata[0] = compspec->data[0];
[729a3c0]394  for (i=1; i < s->fft_size - 1; i++) {
395    s->specdata[i] = compspec->data[i] +
396      I * compspec->data[compspec->length - i];
[d95ff38]397  }
398  s->specdata[s->fft_size - 1] = compspec->data[s->fft_size - 1];
[237f632]399#else
[729a3c0]400  for (i=0; i < s->fft_size; i++) {
401    s->specdata[i] = compspec->data[i];
[d95ff38]402  }
[aadd27a]403#endif
[d95ff38]404  fftw_execute(s->pbw);
[729a3c0]405  for (i = 0; i < output->length; i++) {
406    output->data[i] = s->out[i]*renorm;
407  }
[986131d]408
409#elif defined HAVE_ACCELERATE // using ACCELERATE
[54dd945]410  // convert from real imag  [ r0, r1, ..., rN, iN-1, .., i2, i1]
411  // to vDSP packed format   [ r0, rN, r1, i1, ..., rN-1, iN-1 ]
412  s->out[0] = compspec->data[0];
413  s->out[1] = compspec->data[s->winsize / 2];
414  for (i = 1; i < s->fft_size / 2; i++) {
415    s->out[2 * i] = compspec->data[i];
416    s->out[2 * i + 1] = compspec->data[s->winsize - i];
417  }
418  // convert to split complex format used in vDSP
[39c4721]419  aubio_vDSP_ctoz((aubio_DSPComplex*)s->out, 2, &s->spec, 1, s->fft_size/2);
[54dd945]420  // compute the FFT
[39c4721]421  aubio_vDSP_fft_zrip(s->fftSetup, &s->spec, 1, s->log2fftsize, FFT_INVERSE);
[54dd945]422  // convert result to real output
[39c4721]423  aubio_vDSP_ztoc(&s->spec, 1, (aubio_DSPComplex*)output->data, 2, s->fft_size/2);
[54dd945]424  // apply scaling
425  smpl_t scale = 1.0 / s->winsize;
[39c4721]426  aubio_vDSP_vsmul(output->data, 1, &scale, output->data, 1, s->fft_size);
[986131d]427
428#elif defined HAVE_INTEL_IPP  // using Intel IPP
429
430  // convert from real imag  [ r0, 0, ..., rN, iN-1, .., i2, i1, iN-1] to complex format
431  s->complexOut[0].re = compspec->data[0];
432  s->complexOut[0].im = 0;
433  s->complexOut[s->fft_size / 2].re = compspec->data[s->fft_size / 2];
434  s->complexOut[s->fft_size / 2].im = 0.0;
435  for (i = 1; i < s->fft_size / 2; i++) {
436    s->complexOut[i].re = compspec->data[i];
437    s->complexOut[i].im = compspec->data[s->fft_size - i];
438  }
439  // apply fft
[0ad2e17]440  aubio_ippsFFTInv_CCSToR((const aubio_IppFloat *)s->complexOut, output->data, s->fftSpec, s->memBuffer);
[986131d]441  // apply scaling
[0ad2e17]442  aubio_ippsMulC(output->data, 1.0 / s->winsize, output->data, s->fft_size);
[986131d]443
[54dd945]444#else                         // using OOURA
[7100895]445  smpl_t scale = 2.0 / s->winsize;
[729a3c0]446  s->out[0] = compspec->data[0];
447  s->out[1] = compspec->data[s->winsize / 2];
448  for (i = 1; i < s->fft_size - 1; i++) {
449    s->out[2 * i] = compspec->data[i];
450    s->out[2 * i + 1] = - compspec->data[s->winsize - i];
[aadd27a]451  }
[f50c9503]452  aubio_ooura_rdft(s->winsize, -1, s->out, s->ip, s->w);
[729a3c0]453  for (i=0; i < s->winsize; i++) {
454    output->data[i] = s->out[i] * scale;
455  }
[986131d]456#endif
[237f632]457}
458
[14e300e]459void aubio_fft_get_spectrum(const fvec_t * compspec, cvec_t * spectrum) {
460  aubio_fft_get_phas(compspec, spectrum);
461  aubio_fft_get_norm(compspec, spectrum);
[237f632]462}
463
[14e300e]464void aubio_fft_get_realimag(const cvec_t * spectrum, fvec_t * compspec) {
465  aubio_fft_get_imag(spectrum, compspec);
466  aubio_fft_get_real(spectrum, compspec);
[237f632]467}
468
[14e300e]469void aubio_fft_get_phas(const fvec_t * compspec, cvec_t * spectrum) {
[729a3c0]470  uint_t i;
[d95ff38]471  if (compspec->data[0] < 0) {
472    spectrum->phas[0] = PI;
473  } else {
474    spectrum->phas[0] = 0.;
475  }
[14e300e]476#if defined(HAVE_INTEL_IPP)
477  // convert from real imag  [ r0, r1, ..., rN, iN-1, ..., i2, i1, i0]
478  //                     to  [ r0, r1, ..., rN, i0, i1, i2, ..., iN-1]
479  for (i = 1; i < spectrum->length / 2; i++) {
480    ELEM_SWAP(compspec->data[compspec->length - i],
481        compspec->data[spectrum->length + i - 1]);
482  }
483  aubio_ippsAtan2(compspec->data + spectrum->length,
484      compspec->data + 1, spectrum->phas + 1, spectrum->length - 1);
485  // revert the imaginary part back again
486  for (i = 1; i < spectrum->length / 2; i++) {
487    ELEM_SWAP(compspec->data[spectrum->length + i - 1],
488        compspec->data[compspec->length - i]);
489  }
490#else
[729a3c0]491  for (i=1; i < spectrum->length - 1; i++) {
492    spectrum->phas[i] = ATAN2(compspec->data[compspec->length-i],
493        compspec->data[i]);
[d95ff38]494  }
[14e300e]495#endif
[d95ff38]496  if (compspec->data[compspec->length/2] < 0) {
497    spectrum->phas[spectrum->length - 1] = PI;
498  } else {
499    spectrum->phas[spectrum->length - 1] = 0.;
[aadd27a]500  }
[f88a326]501}
502
[14e300e]503void aubio_fft_get_norm(const fvec_t * compspec, cvec_t * spectrum) {
[729a3c0]504  uint_t i = 0;
[d95ff38]505  spectrum->norm[0] = ABS(compspec->data[0]);
[729a3c0]506  for (i=1; i < spectrum->length - 1; i++) {
507    spectrum->norm[i] = SQRT(SQR(compspec->data[i])
508        + SQR(compspec->data[compspec->length - i]) );
[a47cd35]509  }
[729a3c0]510  spectrum->norm[spectrum->length-1] =
[d95ff38]511    ABS(compspec->data[compspec->length/2]);
[f88a326]512}
513
[14e300e]514void aubio_fft_get_imag(const cvec_t * spectrum, fvec_t * compspec) {
[729a3c0]515  uint_t i;
516  for (i = 1; i < ( compspec->length + 1 ) / 2 /*- 1 + 1*/; i++) {
517    compspec->data[compspec->length - i] =
518      spectrum->norm[i]*SIN(spectrum->phas[i]);
[a47cd35]519  }
[f88a326]520}
521
[14e300e]522void aubio_fft_get_real(const cvec_t * spectrum, fvec_t * compspec) {
[729a3c0]523  uint_t i;
524  for (i = 0; i < compspec->length / 2 + 1; i++) {
[152ed05]525    compspec->data[i] =
[729a3c0]526      spectrum->norm[i]*COS(spectrum->phas[i]);
[aadd27a]527  }
[f88a326]528}
Note: See TracBrowser for help on using the repository browser.