source: src/spectral/dct.c @ 52c1de9

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

src/spectral/dct.c: ooura supports size > 1

  • Property mode set to 100644
File size: 6.0 KB
Line 
1/*
2  Copyright (C) 2018 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/** \file
22
23  Discrete Cosine Transform
24
25  Functions aubio_dct_do() and aubio_dct_rdo() are equivalent to MATLAB/Octave
26  dct() and idct() functions, as well as scipy.fftpack.dct(x, norm='ortho') and
27  scipy.fftpack.idct(x, norm='ortho')
28
29  \example spectral/test-dct.c
30
31*/
32
33#include "aubio_priv.h"
34#include "fvec.h"
35#include "spectral/dct.h"
36
37// function pointers prototypes
38typedef void (*aubio_dct_do_t)(aubio_dct_t * s, const fvec_t * input, fvec_t * output);
39typedef void (*aubio_dct_rdo_t)(aubio_dct_t * s, const fvec_t * input, fvec_t * output);
40typedef void (*del_aubio_dct_t)(aubio_dct_t * s);
41
42#if defined(HAVE_ACCELERATE)
43typedef struct _aubio_dct_opt_t aubio_dct_accelerate_t;
44extern aubio_dct_accelerate_t * new_aubio_dct_accelerate (uint_t size);
45extern void aubio_dct_accelerate_do(aubio_dct_accelerate_t *s, const fvec_t *input, fvec_t *output);
46extern void aubio_dct_accelerate_rdo(aubio_dct_accelerate_t *s, const fvec_t *input, fvec_t *output);
47extern void del_aubio_dct_accelerate (aubio_dct_accelerate_t *s);
48#elif defined(HAVE_FFTW3)
49typedef struct _aubio_dct_opt_t aubio_dct_fftw_t;
50extern aubio_dct_fftw_t * new_aubio_dct_fftw (uint_t size);
51extern void aubio_dct_fftw_do(aubio_dct_fftw_t *s, const fvec_t *input, fvec_t *output);
52extern void aubio_dct_fftw_rdo(aubio_dct_fftw_t *s, const fvec_t *input, fvec_t *output);
53extern void del_aubio_dct_fftw (aubio_dct_fftw_t *s);
54#elif defined(HAVE_INTEL_IPP)
55typedef struct _aubio_dct_ipp_t aubio_dct_ipp_t;
56extern aubio_dct_ipp_t * new_aubio_dct_ipp (uint_t size);
57extern void aubio_dct_ipp_do(aubio_dct_ipp_t *s, const fvec_t *input, fvec_t *output);
58extern void aubio_dct_ipp_rdo(aubio_dct_ipp_t *s, const fvec_t *input, fvec_t *output);
59extern void del_aubio_dct_ipp (aubio_dct_ipp_t *s);
60#else
61typedef struct _aubio_dct_ooura_t aubio_dct_ooura_t;
62extern aubio_dct_ooura_t * new_aubio_dct_ooura (uint_t size);
63extern void aubio_dct_ooura_do(aubio_dct_ooura_t *s, const fvec_t *input, fvec_t *output);
64extern void aubio_dct_ooura_rdo(aubio_dct_ooura_t *s, const fvec_t *input, fvec_t *output);
65extern void del_aubio_dct_ooura (aubio_dct_ooura_t *s);
66#endif
67
68// plain mode
69typedef struct _aubio_dct_plain_t aubio_dct_plain_t;
70extern aubio_dct_plain_t * new_aubio_dct_plain (uint_t size);
71extern void aubio_dct_plain_do(aubio_dct_plain_t *s, const fvec_t *input, fvec_t *output);
72extern void aubio_dct_plain_rdo(aubio_dct_plain_t *s, const fvec_t *input, fvec_t *output);
73extern void del_aubio_dct_plain (aubio_dct_plain_t *s);
74
75struct _aubio_dct_t {
76  void *dct;
77  aubio_dct_do_t dct_do;
78  aubio_dct_rdo_t dct_rdo;
79  del_aubio_dct_t del_dct;
80};
81
82aubio_dct_t* new_aubio_dct (uint_t size) {
83  aubio_dct_t * s = AUBIO_NEW(aubio_dct_t);
84  if ((sint_t)size <= 0) goto beach;
85#if defined(HAVE_ACCELERATE)
86  // TODO check
87  // vDSP supports sizes = f * 2 ** n, where n >= 4 and f in [1, 3, 5, 15]
88  // see https://developer.apple.com/documentation/accelerate/1449930-vdsp_dct_createsetup
89  if (aubio_is_power_of_two(size/16) != 1
90      || (size/16 != 3 && size/16 != 5 && size/16 != 15)) {
91    goto plain;
92  }
93  s->dct = (void *)new_aubio_dct_accelerate (size);
94  if (s->dct) {
95    s->dct_do = (aubio_dct_do_t)aubio_dct_accelerate_do;
96    s->dct_rdo = (aubio_dct_rdo_t)aubio_dct_accelerate_rdo;
97    s->del_dct = (del_aubio_dct_t)del_aubio_dct_accelerate;
98    return s;
99  }
100#elif defined(HAVE_FFTW3)
101  // fftw supports any positive integer size
102  s->dct = (void *)new_aubio_dct_fftw (size);
103  if (s->dct) {
104    s->dct_do = (aubio_dct_do_t)aubio_dct_fftw_do;
105    s->dct_rdo = (aubio_dct_rdo_t)aubio_dct_fftw_rdo;
106    s->del_dct = (del_aubio_dct_t)del_aubio_dct_fftw;
107    return s;
108  } else {
109    AUBIO_WRN("dct: unexcepected error while creating dct_fftw with size %d",
110        size);
111    goto plain;
112  }
113#elif defined(HAVE_INTEL_IPP)
114  // unclear from the docs, but intel ipp seems to support any size
115  s->dct = (void *)new_aubio_dct_ipp (size);
116  if (s->dct) {
117    s->dct_do = (aubio_dct_do_t)aubio_dct_ipp_do;
118    s->dct_rdo = (aubio_dct_rdo_t)aubio_dct_ipp_rdo;
119    s->del_dct = (del_aubio_dct_t)del_aubio_dct_ipp;
120    return s;
121  } else {
122    AUBIO_WRN("dct: unexcepected error while creating dct_ipp with size %d",
123        size);
124    goto plain;
125  }
126#else
127  // ooura support sizes that are power of 2
128  if (aubio_is_power_of_two(size) != 1 || size == 1) {
129    goto plain;
130  }
131  s->dct = (void *)new_aubio_dct_ooura (size);
132  if (s->dct) {
133    s->dct_do = (aubio_dct_do_t)aubio_dct_ooura_do;
134    s->dct_rdo = (aubio_dct_rdo_t)aubio_dct_ooura_rdo;
135    s->del_dct = (del_aubio_dct_t)del_aubio_dct_ooura;
136    return s;
137  }
138#endif
139  // falling back to plain mode
140  AUBIO_WRN("dct: d no optimised implementation could be created for size %d",
141      size);
142plain:
143  s->dct = (void *)new_aubio_dct_plain (size);
144  if (s->dct) {
145    s->dct_do = (aubio_dct_do_t)aubio_dct_plain_do;
146    s->dct_rdo = (aubio_dct_rdo_t)aubio_dct_plain_rdo;
147    s->del_dct = (del_aubio_dct_t)del_aubio_dct_plain;
148    return s;
149  } else {
150    goto beach;
151  }
152beach:
153  AUBIO_ERROR("dct: failed creating with size %d, should be > 0", size);
154  AUBIO_FREE (s);
155  return NULL;
156}
157
158void del_aubio_dct(aubio_dct_t *s) {
159  if (s->dct && s->del_dct) s->del_dct (s->dct);
160  AUBIO_FREE (s);
161}
162
163void aubio_dct_do(aubio_dct_t *s, const fvec_t *input, fvec_t *output) {
164  s->dct_do ((void *)s->dct, input, output);
165}
166
167void aubio_dct_rdo(aubio_dct_t *s, const fvec_t *input, fvec_t *output) {
168  s->dct_rdo ((void *)s->dct, input, output);
169}
Note: See TracBrowser for help on using the repository browser.