source: src/spectral/fft.c @ 1a6ef2c

feature/autosinkfeature/cnnfeature/cnn_orgfeature/constantqfeature/crepefeature/crepe_orgfeature/pitchshiftfeature/pydocstringsfeature/timestretchfix/ffmpeg5pitchshiftsamplertimestretchyinfft+
Last change on this file since 1a6ef2c was d95ff38, checked in by Paul Brossier <piem@piem.org>, 15 years ago

src/spectral: switch to mono

  • Property mode set to 100644
File size: 6.6 KB
Line 
1/*
2  Copyright (C) 2003-2009 Paul Brossier <piem@aubio.org>
3
4  This file is part of aubio.
5
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.
10
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/>.
18
19*/
20
21#include "aubio_priv.h"
22#include "fvec.h"
23#include "cvec.h"
24#include "mathutils.h"
25#include "spectral/fft.h"
26
27/* note that <complex.h> is not included here but only in aubio_priv.h, so that
28 * c++ projects can still use their own complex definition. */
29#include <fftw3.h>
30
31#ifdef HAVE_COMPLEX_H
32#if HAVE_FFTW3F
33/** fft data type with complex.h and fftw3f */
34#define FFTW_TYPE fftwf_complex
35#else
36/** fft data type with complex.h and fftw3 */
37#define FFTW_TYPE fftw_complex
38#endif
39#else
40#if HAVE_FFTW3F
41/** fft data type without complex.h and with fftw3f */
42#define FFTW_TYPE float
43#else
44/** fft data type without complex.h and with fftw */
45#define FFTW_TYPE double
46#endif
47#endif
48
49/** fft data type */
50typedef FFTW_TYPE fft_data_t;
51
52#if HAVE_FFTW3F
53#define fftw_malloc            fftwf_malloc
54#define fftw_free              fftwf_free
55#define fftw_execute           fftwf_execute
56#define fftw_plan_dft_r2c_1d   fftwf_plan_dft_r2c_1d
57#define fftw_plan_dft_c2r_1d   fftwf_plan_dft_c2r_1d
58#define fftw_plan_r2r_1d       fftwf_plan_r2r_1d
59#define fftw_plan              fftwf_plan
60#define fftw_destroy_plan      fftwf_destroy_plan
61#endif
62
63#if HAVE_FFTW3F
64#if HAVE_AUBIO_DOUBLE
65#warning "Using aubio in double precision with fftw3 in single precision"
66#endif /* HAVE_AUBIO_DOUBLE */
67#define real_t float
68#else /* HAVE_FFTW3F */
69#if !HAVE_AUBIO_DOUBLE
70#warning "Using aubio in single precision with fftw3 in double precision"
71#endif /* HAVE_AUBIO_DOUBLE */
72#define real_t double
73#endif /* HAVE_FFTW3F */
74
75struct _aubio_fft_t {
76  uint_t winsize;
77  uint_t fft_size;
78  real_t *in, *out;
79  fftw_plan pfw, pbw;
80  fft_data_t * specdata;     /* complex spectral data */
81  fvec_t * compspec;
82};
83
84aubio_fft_t * new_aubio_fft(uint_t winsize) {
85  aubio_fft_t * s = AUBIO_NEW(aubio_fft_t);
86  uint_t i;
87  s->winsize  = winsize;
88  /* allocate memory */
89  s->in       = AUBIO_ARRAY(real_t,winsize);
90  s->out      = AUBIO_ARRAY(real_t,winsize);
91  s->compspec = new_fvec(winsize);
92  /* create plans */
93#ifdef HAVE_COMPLEX_H
94  s->fft_size = winsize/2 + 1;
95  s->specdata = (fft_data_t*)fftw_malloc(sizeof(fft_data_t)*s->fft_size);
96  s->pfw = fftw_plan_dft_r2c_1d(winsize, s->in,  s->specdata, FFTW_ESTIMATE);
97  s->pbw = fftw_plan_dft_c2r_1d(winsize, s->specdata, s->out, FFTW_ESTIMATE);
98#else
99  s->fft_size = winsize;
100  s->specdata = (fft_data_t*)fftw_malloc(sizeof(fft_data_t)*s->fft_size);
101  s->pfw = fftw_plan_r2r_1d(winsize, s->in,  s->specdata, FFTW_R2HC, FFTW_ESTIMATE);
102  s->pbw = fftw_plan_r2r_1d(winsize, s->specdata, s->out, FFTW_HC2R, FFTW_ESTIMATE);
103#endif
104  for (i = 0; i < s->winsize; i++) {
105    s->in[i] = 0.;
106    s->out[i] = 0.;
107  }
108  for (i = 0; i < s->fft_size; i++) {
109    s->specdata[i] = 0.;
110  }
111  return s;
112}
113
114void del_aubio_fft(aubio_fft_t * s) {
115  /* destroy data */
116  del_fvec(s->compspec);
117  fftw_destroy_plan(s->pfw);
118  fftw_destroy_plan(s->pbw);
119  fftw_free(s->specdata);
120  AUBIO_FREE(s->out);
121  AUBIO_FREE(s->in );
122  AUBIO_FREE(s);
123}
124
125void aubio_fft_do(aubio_fft_t * s, fvec_t * input, cvec_t * spectrum) {
126  aubio_fft_do_complex(s, input, s->compspec);
127  aubio_fft_get_spectrum(s->compspec, spectrum);
128}
129
130void aubio_fft_rdo(aubio_fft_t * s, cvec_t * spectrum, fvec_t * output) {
131  aubio_fft_get_realimag(spectrum, s->compspec);
132  aubio_fft_rdo_complex(s, s->compspec, output);
133}
134
135void aubio_fft_do_complex(aubio_fft_t * s, fvec_t * input, fvec_t * compspec) {
136  uint_t j;
137  for (j=0; j < s->winsize; j++) {
138    s->in[j] = input->data[j];
139  }
140  fftw_execute(s->pfw);
141#ifdef HAVE_COMPLEX_H
142  compspec->data[0] = REAL(s->specdata[0]);
143  for (j = 1; j < s->fft_size -1 ; j++) {
144    compspec->data[j] = REAL(s->specdata[j]);
145    compspec->data[compspec->length - j] = IMAG(s->specdata[j]);
146  }
147  compspec->data[s->fft_size-1] = REAL(s->specdata[s->fft_size-1]);
148#else
149  for (j = 0; j < s->fft_size; j++) {
150    compspec->data[j] = s->specdata[j];
151  }
152#endif
153}
154
155void aubio_fft_rdo_complex(aubio_fft_t * s, fvec_t * compspec, fvec_t * output) {
156  uint_t j;
157  const smpl_t renorm = 1./(smpl_t)s->winsize;
158#ifdef HAVE_COMPLEX_H
159  s->specdata[0] = compspec->data[0];
160  for (j=1; j < s->fft_size - 1; j++) {
161    s->specdata[j] = compspec->data[j] + 
162      I * compspec->data[compspec->length - j];
163  }
164  s->specdata[s->fft_size - 1] = compspec->data[s->fft_size - 1];
165#else
166  for (j=0; j < s->fft_size; j++) {
167    s->specdata[j] = compspec->data[j];
168  }
169#endif
170  fftw_execute(s->pbw);
171  for (j = 0; j < output->length; j++) {
172    output->data[j] = s->out[j]*renorm;
173  }
174}
175
176void aubio_fft_get_spectrum(fvec_t * compspec, cvec_t * spectrum) {
177  aubio_fft_get_phas(compspec, spectrum);
178  aubio_fft_get_norm(compspec, spectrum);
179}
180
181void aubio_fft_get_realimag(cvec_t * spectrum, fvec_t * compspec) {
182  aubio_fft_get_imag(spectrum, compspec);
183  aubio_fft_get_real(spectrum, compspec);
184}
185
186void aubio_fft_get_phas(fvec_t * compspec, cvec_t * spectrum) {
187  uint_t j;
188  if (compspec->data[0] < 0) {
189    spectrum->phas[0] = PI;
190  } else {
191    spectrum->phas[0] = 0.;
192  }
193  for (j=1; j < spectrum->length - 1; j++) {
194    spectrum->phas[j] = ATAN2(compspec->data[compspec->length-j],
195        compspec->data[j]);
196  }
197  if (compspec->data[compspec->length/2] < 0) {
198    spectrum->phas[spectrum->length - 1] = PI;
199  } else {
200    spectrum->phas[spectrum->length - 1] = 0.;
201  }
202}
203
204void aubio_fft_get_norm(fvec_t * compspec, cvec_t * spectrum) {
205  uint_t j = 0;
206  spectrum->norm[0] = ABS(compspec->data[0]);
207  for (j=1; j < spectrum->length - 1; j++) {
208    spectrum->norm[j] = SQRT(SQR(compspec->data[j]) 
209        + SQR(compspec->data[compspec->length - j]) );
210  }
211  spectrum->norm[spectrum->length-1] = 
212    ABS(compspec->data[compspec->length/2]);
213}
214
215void aubio_fft_get_imag(cvec_t * spectrum, fvec_t * compspec) {
216  uint_t j;
217  for (j = 1; j < ( compspec->length + 1 ) / 2 /*- 1 + 1*/; j++) {
218    compspec->data[compspec->length - j] =
219      spectrum->norm[j]*SIN(spectrum->phas[j]);
220  }
221}
222
223void aubio_fft_get_real(cvec_t * spectrum, fvec_t * compspec) {
224  uint_t j;
225  for (j = 0; j < compspec->length / 2 + 1; j++) {
226    compspec->data[j] = 
227      spectrum->norm[j]*COS(spectrum->phas[j]);
228  }
229}
Note: See TracBrowser for help on using the repository browser.