matlabfunctions.cpp 9.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321
  1. //-----------------------------------------------------------------------------
  2. // Copyright 2012 Masanori Morise
  3. // Author: mmorise [at] yamanashi.ac.jp (Masanori Morise)
  4. // Last update: 2017/02/01
  5. //
  6. // Matlab functions implemented for WORLD
  7. // Since these functions are implemented as the same function of Matlab,
  8. // the source code does not follow the style guide (Names of variables
  9. // and functions).
  10. // Please see the reference of Matlab to show the usage of functions.
  11. // Caution:
  12. // The functions wavread() and wavwrite() were removed to the /src.
  13. // they were moved to the test/audioio.cpp. (2016/01/28)
  14. //-----------------------------------------------------------------------------
  15. #include "world/matlabfunctions.h"
  16. #include <math.h>
  17. #include <stdint.h>
  18. #include "world/constantnumbers.h"
  19. namespace {
  20. //-----------------------------------------------------------------------------
  21. // FilterForDecimate() calculates the coefficients of low-pass filter and
  22. // carries out the filtering. This function is only used for decimate().
  23. //-----------------------------------------------------------------------------
  24. static void FilterForDecimate(const double *x, int x_length, int r, double *y) {
  25. double a[3], b[2]; // filter Coefficients
  26. switch (r) {
  27. case 11: // fs : 44100 (default)
  28. a[0] = 2.450743295230728;
  29. a[1] = -2.06794904601978;
  30. a[2] = 0.59574774438332101;
  31. b[0] = 0.0026822508007163792;
  32. b[1] = 0.0080467524021491377;
  33. break;
  34. case 12: // fs : 48000
  35. a[0] = 2.4981398605924205;
  36. a[1] = -2.1368928194784025;
  37. a[2] = 0.62187513816221485;
  38. b[0] = 0.0021097275904709001;
  39. b[1] = 0.0063291827714127002;
  40. break;
  41. case 10:
  42. a[0] = 2.3936475118069387;
  43. a[1] = -1.9873904075111861;
  44. a[2] = 0.5658879979027055;
  45. b[0] = 0.0034818622251927556;
  46. b[1] = 0.010445586675578267;
  47. break;
  48. case 9:
  49. a[0] = 2.3236003491759578;
  50. a[1] = -1.8921545617463598;
  51. a[2] = 0.53148928133729068;
  52. b[0] = 0.0046331164041389372;
  53. b[1] = 0.013899349212416812;
  54. break;
  55. case 8: // fs : 32000
  56. a[0] = 2.2357462340187593;
  57. a[1] = -1.7780899984041358;
  58. a[2] = 0.49152555365968692;
  59. b[0] = 0.0063522763407111993;
  60. b[1] = 0.019056829022133598;
  61. break;
  62. case 7:
  63. a[0] = 2.1225239019534703;
  64. a[1] = -1.6395144861046302;
  65. a[2] = 0.44469707800587366;
  66. b[0] = 0.0090366882681608418;
  67. b[1] = 0.027110064804482525;
  68. break;
  69. case 6: // fs : 24000 and 22050
  70. a[0] = 1.9715352749512141;
  71. a[1] = -1.4686795689225347;
  72. a[2] = 0.3893908434965701;
  73. b[0] = 0.013469181309343825;
  74. b[1] = 0.040407543928031475;
  75. break;
  76. case 5:
  77. a[0] = 1.7610939654280557;
  78. a[1] = -1.2554914843859768;
  79. a[2] = 0.3237186507788215;
  80. b[0] = 0.021334858522387423;
  81. b[1] = 0.06400457556716227;
  82. break;
  83. case 4: // fs : 16000
  84. a[0] = 1.4499664446880227;
  85. a[1] = -0.98943497080950582;
  86. a[2] = 0.24578252340690215;
  87. b[0] = 0.036710750339322612;
  88. b[1] = 0.11013225101796784;
  89. break;
  90. case 3:
  91. a[0] = 0.95039378983237421;
  92. a[1] = -0.67429146741526791;
  93. a[2] = 0.15412211621346475;
  94. b[0] = 0.071221945171178636;
  95. b[1] = 0.21366583551353591;
  96. break;
  97. case 2: // fs : 8000
  98. a[0] = 0.041156734567757189;
  99. a[1] = -0.42599112459189636;
  100. a[2] = 0.041037215479961225;
  101. b[0] = 0.16797464681802227;
  102. b[1] = 0.50392394045406674;
  103. break;
  104. default:
  105. a[0] = 0.0;
  106. a[1] = 0.0;
  107. a[2] = 0.0;
  108. b[0] = 0.0;
  109. b[1] = 0.0;
  110. }
  111. // Filtering on time domain.
  112. double w[3] = {0.0, 0.0, 0.0};
  113. double wt;
  114. for (int i = 0; i < x_length; ++i) {
  115. wt = x[i] + a[0] * w[0] + a[1] * w[1] + a[2] * w[2];
  116. y[i] = b[0] * wt + b[1] * w[0] + b[1] * w[1] + b[0] * w[2];
  117. w[2] = w[1];
  118. w[1] = w[0];
  119. w[0] = wt;
  120. }
  121. }
  122. } // namespace
  123. void fftshift(const double *x, int x_length, double *y) {
  124. for (int i = 0; i < x_length / 2; ++i) {
  125. y[i] = x[i + x_length / 2];
  126. y[i + x_length / 2] = x[i];
  127. }
  128. }
  129. void histc(const double *x, int x_length, const double *edges,
  130. int edges_length, int *index) {
  131. int count = 1;
  132. int i = 0;
  133. for (; i < edges_length; ++i) {
  134. index[i] = 1;
  135. if (edges[i] >= x[0]) break;
  136. }
  137. for (; i < edges_length; ++i) {
  138. if (edges[i] < x[count]) {
  139. index[i] = count;
  140. } else {
  141. index[i--] = count++;
  142. }
  143. if (count == x_length) break;
  144. }
  145. count--;
  146. for (i++; i < edges_length; ++i) index[i] = count;
  147. }
  148. void interp1(const double *x, const double *y, int x_length, const double *xi,
  149. int xi_length, double *yi) {
  150. double *h = new double[x_length - 1];
  151. int *k = new int[xi_length];
  152. for (int i = 0; i < x_length - 1; ++i) h[i] = x[i + 1] - x[i];
  153. for (int i = 0; i < xi_length; ++i) {
  154. k[i] = 0;
  155. }
  156. histc(x, x_length, xi, xi_length, k);
  157. for (int i = 0; i < xi_length; ++i) {
  158. double s = (xi[i] - x[k[i] - 1]) / h[k[i] - 1];
  159. yi[i] = y[k[i] - 1] + s * (y[k[i]] - y[k[i] - 1]);
  160. }
  161. delete[] k;
  162. delete[] h;
  163. }
  164. void decimate(const double *x, int x_length, int r, double *y) {
  165. const int kNFact = 9;
  166. double *tmp1 = new double[x_length + kNFact * 2];
  167. double *tmp2 = new double[x_length + kNFact * 2];
  168. for (int i = 0; i < kNFact; ++i) tmp1[i] = 2 * x[0] - x[kNFact - i];
  169. for (int i = kNFact; i < kNFact + x_length; ++i) tmp1[i] = x[i - kNFact];
  170. for (int i = kNFact + x_length; i < 2 * kNFact + x_length; ++i)
  171. tmp1[i] = 2 * x[x_length - 1] - x[x_length - 2 - (i - (kNFact + x_length))];
  172. FilterForDecimate(tmp1, 2 * kNFact + x_length, r, tmp2);
  173. for (int i = 0; i < 2 * kNFact + x_length; ++i)
  174. tmp1[i] = tmp2[2 * kNFact + x_length - i - 1];
  175. FilterForDecimate(tmp1, 2 * kNFact + x_length, r, tmp2);
  176. for (int i = 0; i < 2 * kNFact + x_length; ++i)
  177. tmp1[i] = tmp2[2 * kNFact + x_length - i - 1];
  178. int nout = (x_length - 1) / r + 1;
  179. int nbeg = r - r * nout + x_length;
  180. int count = 0;
  181. for (int i = nbeg; i < x_length + kNFact; i += r)
  182. y[count++] = tmp1[i + kNFact - 1];
  183. delete[] tmp1;
  184. delete[] tmp2;
  185. }
  186. int matlab_round(double x) {
  187. return x > 0 ? static_cast<int>(x + 0.5) : static_cast<int>(x - 0.5);
  188. }
  189. void diff(const double *x, int x_length, double *y) {
  190. for (int i = 0; i < x_length - 1; ++i) y[i] = x[i + 1] - x[i];
  191. }
  192. void interp1Q(double x, double shift, const double *y, int x_length,
  193. const double *xi, int xi_length, double *yi) {
  194. double *xi_fraction = new double[xi_length];
  195. double *delta_y = new double[x_length];
  196. int *xi_base = new int[xi_length];
  197. double delta_x = shift;
  198. for (int i = 0; i < xi_length; ++i) {
  199. xi_base[i] = static_cast<int>((xi[i] - x) / delta_x);
  200. xi_fraction[i] = (xi[i] - x) / delta_x - xi_base[i];
  201. }
  202. diff(y, x_length, delta_y);
  203. delta_y[x_length - 1] = 0.0;
  204. for (int i = 0; i < xi_length; ++i)
  205. yi[i] = y[xi_base[i]] + delta_y[xi_base[i]] * xi_fraction[i];
  206. // Bug was fixed at 2013/07/14 by M. Morise
  207. delete[] xi_fraction;
  208. delete[] xi_base;
  209. delete[] delta_y;
  210. }
  211. // You must not use these variables.
  212. // Note:
  213. // I have no idea to implement the randn() and randn_reseed() without the
  214. // global variables. If you have a good idea, please give me the information.
  215. static uint32_t g_randn_x = 123456789;
  216. static uint32_t g_randn_y = 362436069;
  217. static uint32_t g_randn_z = 521288629;
  218. static uint32_t g_randn_w = 88675123;
  219. void randn_reseed() {
  220. g_randn_x = 123456789;
  221. g_randn_y = 362436069;
  222. g_randn_z = 521288629;
  223. g_randn_w = 88675123;
  224. }
  225. double randn(void) {
  226. uint32_t t;
  227. t = g_randn_x ^ (g_randn_x << 11);
  228. g_randn_x = g_randn_y;
  229. g_randn_y = g_randn_z;
  230. g_randn_z = g_randn_w;
  231. g_randn_w = (g_randn_w ^ (g_randn_w >> 19)) ^ (t ^ (t >> 8));
  232. uint32_t tmp = g_randn_w >> 4;
  233. for (int i = 0; i < 11; ++i) {
  234. t = g_randn_x ^ (g_randn_x << 11);
  235. g_randn_x = g_randn_y;
  236. g_randn_y = g_randn_z;
  237. g_randn_z = g_randn_w;
  238. g_randn_w = (g_randn_w ^ (g_randn_w >> 19)) ^ (t ^ (t >> 8));
  239. tmp += g_randn_w >> 4;
  240. }
  241. return tmp / 268435456.0 - 6.0;
  242. }
  243. void fast_fftfilt(const double *x, int x_length, const double *h, int h_length,
  244. int fft_size, const ForwardRealFFT *forward_real_fft,
  245. const InverseRealFFT *inverse_real_fft, double *y) {
  246. fft_complex *x_spectrum = new fft_complex[fft_size];
  247. for (int i = 0; i < x_length; ++i)
  248. forward_real_fft->waveform[i] = x[i] / fft_size;
  249. for (int i = x_length; i < fft_size; ++i)
  250. forward_real_fft->waveform[i] = 0.0;
  251. fft_execute(forward_real_fft->forward_fft);
  252. for (int i = 0; i <= fft_size / 2; ++i) {
  253. x_spectrum[i][0] = forward_real_fft->spectrum[i][0];
  254. x_spectrum[i][1] = forward_real_fft->spectrum[i][1];
  255. }
  256. for (int i = 0; i < h_length; ++i)
  257. forward_real_fft->waveform[i] = h[i] / fft_size;
  258. for (int i = h_length; i < fft_size; ++i)
  259. forward_real_fft->waveform[i] = 0.0;
  260. fft_execute(forward_real_fft->forward_fft);
  261. for (int i = 0; i <= fft_size / 2; ++i) {
  262. inverse_real_fft->spectrum[i][0] =
  263. x_spectrum[i][0] * forward_real_fft->spectrum[i][0] -
  264. x_spectrum[i][1] * forward_real_fft->spectrum[i][1];
  265. inverse_real_fft->spectrum[i][1] =
  266. x_spectrum[i][0] * forward_real_fft->spectrum[i][1] +
  267. x_spectrum[i][1] * forward_real_fft->spectrum[i][0];
  268. }
  269. fft_execute(inverse_real_fft->inverse_fft);
  270. for (int i = 0; i < fft_size; ++i)
  271. y[i] = inverse_real_fft->waveform[i];
  272. delete[] x_spectrum;
  273. }
  274. double matlab_std(const double *x, int x_length) {
  275. double average = 0.0;
  276. for (int i = 0; i < x_length; ++i) average += x[i];
  277. average /= x_length;
  278. double s = 0.0;
  279. for (int i = 0; i < x_length; ++i) s += pow(x[i] - average, 2.0);
  280. s /= (x_length - 1);
  281. return sqrt(s);
  282. }