annotate hmm/hmm.c @ 30:a251fb0de594

* Make MFCC able to accept already-FFT'd input, and simplify API a bit * Add log power value to MFCC, restore windowing, and avoid some heap allocs * In HMM, bail out of iteration if loglik hits NaN
author cannam
date Fri, 18 Jan 2008 13:24:12 +0000
parents d3a856b44c43
children 00603b8a940f
rev   line source
cannam@19 1 /*
cannam@19 2 * hmm.c
cannam@19 3 * soundbite
cannam@19 4 *
cannam@19 5 * Created by Mark Levy on 12/02/2006.
cannam@19 6 * Copyright 2006 Centre for Digital Music, Queen Mary, University of London. All rights reserved.
cannam@19 7 *
cannam@19 8 */
cannam@19 9
cannam@19 10 #include <stdio.h>
cannam@19 11 #include <math.h>
cannam@19 12 #include <stdlib.h>
cannam@19 13 #include <float.h>
cannam@19 14 #include <time.h> /* to seed random number generator */
cannam@19 15 #include <clapack.h> /* LAPACK for matrix inversion */
cannam@19 16 #ifdef _MAC_OS_X
cannam@19 17 #include <vecLib/cblas.h>
cannam@19 18 #else
cannam@19 19 #include <cblas.h> /* BLAS for matrix multiplication */
cannam@19 20 #endif
cannam@19 21
cannam@19 22 #include "hmm.h"
cannam@19 23
cannam@19 24 model_t* hmm_init(double** x, int T, int L, int N)
cannam@19 25 {
cannam@19 26 int i, j, d, e, t;
cannam@19 27 double s, ss;
cannam@19 28
cannam@19 29 model_t* model;
cannam@19 30 model = (model_t*) malloc(sizeof(model_t));
cannam@19 31 model->N = N;
cannam@19 32 model->L = L;
cannam@19 33 model->p0 = (double*) malloc(N*sizeof(double));
cannam@19 34 model->a = (double**) malloc(N*sizeof(double*));
cannam@19 35 model->mu = (double**) malloc(N*sizeof(double*));
cannam@19 36 for (i = 0; i < N; i++)
cannam@19 37 {
cannam@19 38 model->a[i] = (double*) malloc(N*sizeof(double));
cannam@19 39 model->mu[i] = (double*) malloc(L*sizeof(double));
cannam@19 40 }
cannam@19 41 model->cov = (double**) malloc(L*sizeof(double*));
cannam@19 42 for (i = 0; i < L; i++)
cannam@19 43 model->cov[i] = (double*) malloc(L*sizeof(double));
cannam@19 44
cannam@19 45 srand(time(0));
cannam@19 46 double* global_mean = (double*) malloc(L*sizeof(double));
cannam@19 47
cannam@19 48 /* find global mean */
cannam@19 49 for (d = 0; d < L; d++)
cannam@19 50 {
cannam@19 51 global_mean[d] = 0;
cannam@19 52 for (t = 0; t < T; t++)
cannam@19 53 global_mean[d] += x[t][d];
cannam@19 54 global_mean[d] /= T;
cannam@19 55 }
cannam@19 56
cannam@19 57 /* calculate global diagonal covariance */
cannam@19 58 for (d = 0; d < L; d++)
cannam@19 59 {
cannam@19 60 for (e = 0; e < L; e++)
cannam@19 61 model->cov[d][e] = 0;
cannam@19 62 for (t = 0; t < T; t++)
cannam@19 63 model->cov[d][d] += (x[t][d] - global_mean[d]) * (x[t][d] - global_mean[d]);
cannam@19 64 model->cov[d][d] /= T-1;
cannam@19 65 }
cannam@19 66
cannam@19 67 /* set all means close to global mean */
cannam@19 68 for (i = 0; i < N; i++)
cannam@19 69 {
cannam@19 70 for (d = 0; d < L; d++)
cannam@19 71 {
cannam@19 72 /* add some random noise related to covariance */
cannam@19 73 /* ideally the random number would be Gaussian(0,1), as a hack we make it uniform on [-0.25,0.25] */
cannam@19 74 model->mu[i][d] = global_mean[d] + (0.5 * rand() / (double) RAND_MAX - 0.25) * sqrt(model->cov[d][d]);
cannam@19 75 }
cannam@19 76 }
cannam@19 77
cannam@19 78 /* random intial and transition probs */
cannam@19 79 s = 0;
cannam@19 80 for (i = 0; i < N; i++)
cannam@19 81 {
cannam@19 82 model->p0[i] = 1 + rand() / (double) RAND_MAX;
cannam@19 83 s += model->p0[i];
cannam@19 84 ss = 0;
cannam@19 85 for (j = 0; j < N; j++)
cannam@19 86 {
cannam@19 87 model->a[i][j] = 1 + rand() / (double) RAND_MAX;
cannam@19 88 ss += model->a[i][j];
cannam@19 89 }
cannam@19 90 for (j = 0; j < N; j++)
cannam@19 91 {
cannam@19 92 model->a[i][j] /= ss;
cannam@19 93 }
cannam@19 94 }
cannam@19 95 for (i = 0; i < N; i++)
cannam@19 96 model->p0[i] /= s;
cannam@19 97
cannam@19 98 free(global_mean);
cannam@19 99
cannam@19 100 return model;
cannam@19 101 }
cannam@19 102
cannam@19 103 void hmm_close(model_t* model)
cannam@19 104 {
cannam@19 105 int i;
cannam@19 106
cannam@19 107 for (i = 0; i < model->N; i++)
cannam@19 108 {
cannam@19 109 free(model->a[i]);
cannam@19 110 free(model->mu[i]);
cannam@19 111 }
cannam@19 112 free(model->a);
cannam@19 113 free(model->mu);
cannam@19 114 for (i = 0; i < model->L; i++)
cannam@19 115 free(model->cov[i]);
cannam@19 116 free(model->cov);
cannam@19 117 free(model);
cannam@19 118 }
cannam@19 119
cannam@19 120 void hmm_train(double** x, int T, model_t* model)
cannam@19 121 {
cannam@19 122 int i, t;
cannam@19 123 double loglik; /* overall log-likelihood at each iteration */
cannam@19 124
cannam@19 125 int N = model->N;
cannam@19 126 int L = model->L;
cannam@19 127 double* p0 = model->p0;
cannam@19 128 double** a = model->a;
cannam@19 129 double** mu = model->mu;
cannam@19 130 double** cov = model->cov;
cannam@19 131
cannam@19 132 /* allocate memory */
cannam@19 133 double** gamma = (double**) malloc(T*sizeof(double*));
cannam@19 134 double*** xi = (double***) malloc(T*sizeof(double**));
cannam@19 135 for (t = 0; t < T; t++)
cannam@19 136 {
cannam@19 137 gamma[t] = (double*) malloc(N*sizeof(double));
cannam@19 138 xi[t] = (double**) malloc(N*sizeof(double*));
cannam@19 139 for (i = 0; i < N; i++)
cannam@19 140 xi[t][i] = (double*) malloc(N*sizeof(double));
cannam@19 141 }
cannam@19 142
cannam@19 143 /* temporary memory */
cannam@19 144 double* gauss_y = (double*) malloc(L*sizeof(double));
cannam@19 145 double* gauss_z = (double*) malloc(L*sizeof(double));
cannam@19 146
cannam@19 147 /* obs probs P(j|{x}) */
cannam@19 148 double** b = (double**) malloc(T*sizeof(double*));
cannam@19 149 for (t = 0; t < T; t++)
cannam@19 150 b[t] = (double*) malloc(N*sizeof(double));
cannam@19 151
cannam@19 152 /* inverse covariance and its determinant */
cannam@19 153 double** icov = (double**) malloc(L*sizeof(double*));
cannam@19 154 for (i = 0; i < L; i++)
cannam@19 155 icov[i] = (double*) malloc(L*sizeof(double));
cannam@19 156 double detcov;
cannam@19 157
cannam@19 158 double thresh = 0.0001;
cannam@19 159 int niter = 50;
cannam@19 160 int iter = 0;
cannam@19 161 double loglik1, loglik2;
cannam@30 162 int foundnan = 0;
cannam@30 163
cannam@30 164 while (iter < niter && !foundnan && !(iter > 1 && (loglik - loglik1) < thresh * (loglik1 - loglik2)))
cannam@19 165 {
cannam@19 166 ++iter;
cannam@19 167
cannam@19 168 fprintf(stderr, "calculating obsprobs...\n");
cannam@19 169 fflush(stderr);
cannam@19 170
cannam@19 171 /* precalculate obs probs */
cannam@19 172 invert(cov, L, icov, &detcov);
cannam@19 173
cannam@19 174 for (t = 0; t < T; t++)
cannam@19 175 {
cannam@19 176 //int allzero = 1;
cannam@19 177 for (i = 0; i < N; i++)
cannam@19 178 {
cannam@19 179 b[t][i] = exp(loggauss(x[t], L, mu[i], icov, detcov, gauss_y, gauss_z));
cannam@19 180
cannam@19 181 //if (b[t][i] != 0)
cannam@19 182 // allzero = 0;
cannam@19 183 }
cannam@19 184 /*
cannam@19 185 if (allzero)
cannam@19 186 {
cannam@19 187 printf("all the b[t][i] were zero for t = %d, correcting...\n", t);
cannam@19 188 for (i = 0; i < N; i++)
cannam@19 189 {
cannam@19 190 b[t][i] = 0.00001;
cannam@19 191 }
cannam@19 192 }
cannam@19 193 */
cannam@19 194 }
cannam@19 195
cannam@19 196 fprintf(stderr, "forwards-backwards...\n");
cannam@19 197 fflush(stderr);
cannam@19 198
cannam@19 199 forward_backwards(xi, gamma, &loglik, &loglik1, &loglik2, iter, N, T, p0, a, b);
cannam@19 200
cannam@19 201 fprintf(stderr, "iteration %d: loglik = %f\n", iter, loglik);
cannam@19 202 fprintf(stderr, "re-estimation...\n");
cannam@19 203 fflush(stderr);
cannam@30 204
cannam@30 205 if (isnan(loglik)) {
cannam@30 206 foundnan = 1;
cannam@30 207 continue;
cannam@30 208 }
cannam@19 209
cannam@19 210 baum_welch(p0, a, mu, cov, N, T, L, x, xi, gamma);
cannam@19 211
cannam@19 212 /*
cannam@19 213 printf("a:\n");
cannam@19 214 for (i = 0; i < model->N; i++)
cannam@19 215 {
cannam@19 216 for (j = 0; j < model->N; j++)
cannam@19 217 printf("%f ", model->a[i][j]);
cannam@19 218 printf("\n");
cannam@19 219 }
cannam@19 220 printf("\n\n");
cannam@19 221 */
cannam@19 222 //hmm_print(model);
cannam@19 223 }
cannam@19 224
cannam@19 225 /* deallocate memory */
cannam@19 226 for (t = 0; t < T; t++)
cannam@19 227 {
cannam@19 228 free(gamma[t]);
cannam@19 229 free(b[t]);
cannam@19 230 for (i = 0; i < N; i++)
cannam@19 231 free(xi[t][i]);
cannam@19 232 free(xi[t]);
cannam@19 233 }
cannam@19 234 free(gamma);
cannam@19 235 free(xi);
cannam@19 236 free(b);
cannam@19 237
cannam@19 238 for (i = 0; i < L; i++)
cannam@19 239 free(icov[i]);
cannam@19 240 free(icov);
cannam@19 241
cannam@19 242 free(gauss_y);
cannam@19 243 free(gauss_z);
cannam@19 244 }
cannam@19 245
cannam@19 246 void mlss_reestimate(double* p0, double** a, double** mu, double** cov, int N, int T, int L, int* q, double** x)
cannam@19 247 {
cannam@19 248 /* fit a single Gaussian to observations in each state */
cannam@19 249
cannam@19 250 /* calculate the mean observation in each state */
cannam@19 251
cannam@19 252 /* calculate the overall covariance */
cannam@19 253
cannam@19 254 /* count transitions */
cannam@19 255
cannam@19 256 /* estimate initial probs from transitions (???) */
cannam@19 257 }
cannam@19 258
cannam@19 259 void baum_welch(double* p0, double** a, double** mu, double** cov, int N, int T, int L, double** x, double*** xi, double** gamma)
cannam@19 260 {
cannam@19 261 int i, j, t;
cannam@19 262
cannam@19 263 double* sum_gamma = (double*) malloc(N*sizeof(double));
cannam@19 264
cannam@19 265 /* temporary memory */
cannam@19 266 double* u = (double*) malloc(L*L*sizeof(double));
cannam@19 267 double* yy = (double*) malloc(T*L*sizeof(double));
cannam@19 268 double* yy2 = (double*) malloc(T*L*sizeof(double));
cannam@19 269
cannam@19 270 /* re-estimate transition probs */
cannam@19 271 for (i = 0; i < N; i++)
cannam@19 272 {
cannam@19 273 sum_gamma[i] = 0;
cannam@19 274 for (t = 0; t < T-1; t++)
cannam@19 275 sum_gamma[i] += gamma[t][i];
cannam@19 276 }
cannam@19 277
cannam@19 278 for (i = 0; i < N; i++)
cannam@19 279 {
cannam@19 280 if (sum_gamma[i] == 0)
cannam@19 281 {
cannam@19 282 fprintf(stderr, "sum_gamma[%d] was zero...\n", i);
cannam@19 283 }
cannam@19 284 //double s = 0;
cannam@19 285 for (j = 0; j < N; j++)
cannam@19 286 {
cannam@19 287 a[i][j] = 0;
cannam@30 288 if (sum_gamma[i] == 0.) continue;
cannam@19 289 for (t = 0; t < T-1; t++)
cannam@19 290 a[i][j] += xi[t][i][j];
cannam@19 291 //s += a[i][j];
cannam@19 292 a[i][j] /= sum_gamma[i];
cannam@19 293 }
cannam@19 294 /*
cannam@19 295 for (j = 0; j < N; j++)
cannam@19 296 {
cannam@19 297 a[i][j] /= s;
cannam@19 298 }
cannam@19 299 */
cannam@19 300 }
cannam@19 301
cannam@19 302 /* NB: now we need to sum gamma over all t */
cannam@19 303 for (i = 0; i < N; i++)
cannam@19 304 sum_gamma[i] += gamma[T-1][i];
cannam@19 305
cannam@19 306 /* re-estimate initial probs */
cannam@19 307 for (i = 0; i < N; i++)
cannam@19 308 p0[i] = gamma[0][i];
cannam@19 309
cannam@19 310 /* re-estimate covariance */
cannam@19 311 int d, e;
cannam@19 312 double sum_sum_gamma = 0;
cannam@19 313 for (i = 0; i < N; i++)
cannam@19 314 sum_sum_gamma += sum_gamma[i];
cannam@19 315
cannam@19 316 /*
cannam@19 317 for (d = 0; d < L; d++)
cannam@19 318 {
cannam@19 319 for (e = d; e < L; e++)
cannam@19 320 {
cannam@19 321 cov[d][e] = 0;
cannam@19 322 for (t = 0; t < T; t++)
cannam@19 323 for (j = 0; j < N; j++)
cannam@19 324 cov[d][e] += gamma[t][j] * (x[t][d] - mu[j][d]) * (x[t][e] - mu[j][e]);
cannam@19 325
cannam@19 326 cov[d][e] /= sum_sum_gamma;
cannam@19 327
cannam@19 328 if (isnan(cov[d][e]))
cannam@19 329 {
cannam@19 330 printf("cov[%d][%d] was nan\n", d, e);
cannam@19 331 for (j = 0; j < N; j++)
cannam@19 332 for (i = 0; i < L; i++)
cannam@19 333 if (isnan(mu[j][i]))
cannam@19 334 printf("mu[%d][%d] was nan\n", j, i);
cannam@19 335 for (t = 0; t < T; t++)
cannam@19 336 for (j = 0; j < N; j++)
cannam@19 337 if (isnan(gamma[t][j]))
cannam@19 338 printf("gamma[%d][%d] was nan\n", t, j);
cannam@19 339 exit(-1);
cannam@19 340 }
cannam@19 341 }
cannam@19 342 }
cannam@19 343 for (d = 0; d < L; d++)
cannam@19 344 for (e = 0; e < d; e++)
cannam@19 345 cov[d][e] = cov[e][d];
cannam@19 346 */
cannam@19 347
cannam@19 348 /* using BLAS */
cannam@19 349 for (d = 0; d < L; d++)
cannam@19 350 for (e = 0; e < L; e++)
cannam@19 351 cov[d][e] = 0;
cannam@19 352
cannam@19 353 for (j = 0; j < N; j++)
cannam@19 354 {
cannam@19 355 for (d = 0; d < L; d++)
cannam@19 356 for (t = 0; t < T; t++)
cannam@19 357 {
cannam@19 358 yy[d*T+t] = x[t][d] - mu[j][d];
cannam@19 359 yy2[d*T+t] = gamma[t][j] * (x[t][d] - mu[j][d]);
cannam@19 360 }
cannam@19 361
cannam@19 362 cblas_dgemm(CblasColMajor, CblasTrans, CblasNoTrans, L, L, T, 1.0, yy, T, yy2, T, 0, u, L);
cannam@19 363
cannam@19 364 for (e = 0; e < L; e++)
cannam@19 365 for (d = 0; d < L; d++)
cannam@19 366 cov[d][e] += u[e*L+d];
cannam@19 367 }
cannam@19 368
cannam@19 369 for (d = 0; d < L; d++)
cannam@19 370 for (e = 0; e < L; e++)
cannam@19 371 cov[d][e] /= T; /* sum_sum_gamma; */
cannam@19 372
cannam@19 373 //printf("sum_sum_gamma = %f\n", sum_sum_gamma); /* fine, = T IS THIS ALWAYS TRUE with pooled cov?? */
cannam@19 374
cannam@19 375 /* re-estimate means */
cannam@19 376 for (j = 0; j < N; j++)
cannam@19 377 {
cannam@19 378 for (d = 0; d < L; d++)
cannam@19 379 {
cannam@19 380 mu[j][d] = 0;
cannam@19 381 for (t = 0; t < T; t++)
cannam@19 382 mu[j][d] += gamma[t][j] * x[t][d];
cannam@19 383 mu[j][d] /= sum_gamma[j];
cannam@19 384 }
cannam@19 385 }
cannam@19 386
cannam@19 387 /* deallocate memory */
cannam@19 388 free(sum_gamma);
cannam@19 389 free(yy);
cannam@19 390 free(yy2);
cannam@19 391 free(u);
cannam@19 392 }
cannam@19 393
cannam@19 394 void forward_backwards(double*** xi, double** gamma, double* loglik, double* loglik1, double* loglik2, int iter, int N, int T, double* p0, double** a, double** b)
cannam@19 395 {
cannam@19 396 /* forwards-backwards with scaling */
cannam@19 397 int i, j, t;
cannam@19 398
cannam@19 399 double** alpha = (double**) malloc(T*sizeof(double*));
cannam@19 400 double** beta = (double**) malloc(T*sizeof(double*));
cannam@19 401 for (t = 0; t < T; t++)
cannam@19 402 {
cannam@19 403 alpha[t] = (double*) malloc(N*sizeof(double));
cannam@19 404 beta[t] = (double*) malloc(N*sizeof(double));
cannam@19 405 }
cannam@19 406
cannam@19 407 /* scaling coefficients */
cannam@19 408 double* c = (double*) malloc(T*sizeof(double));
cannam@19 409
cannam@19 410 /* calculate forward probs and scale coefficients */
cannam@19 411 c[0] = 0;
cannam@19 412 for (i = 0; i < N; i++)
cannam@19 413 {
cannam@19 414 alpha[0][i] = p0[i] * b[0][i];
cannam@19 415 c[0] += alpha[0][i];
cannam@19 416
cannam@19 417 //printf("p0[%d] = %f, b[0][%d] = %f\n", i, p0[i], i, b[0][i]);
cannam@19 418 }
cannam@19 419 c[0] = 1 / c[0];
cannam@19 420 for (i = 0; i < N; i++)
cannam@19 421 {
cannam@19 422 alpha[0][i] *= c[0];
cannam@19 423
cannam@19 424 //printf("alpha[0][%d] = %f\n", i, alpha[0][i]); /* OK agrees with Matlab */
cannam@19 425 }
cannam@19 426
cannam@19 427 *loglik1 = *loglik;
cannam@19 428 *loglik = -log(c[0]);
cannam@19 429 if (iter == 2)
cannam@19 430 *loglik2 = *loglik;
cannam@19 431
cannam@19 432 for (t = 1; t < T; t++)
cannam@19 433 {
cannam@19 434 c[t] = 0;
cannam@19 435 for (j = 0; j < N; j++)
cannam@19 436 {
cannam@19 437 alpha[t][j] = 0;
cannam@19 438 for (i = 0; i < N; i++)
cannam@19 439 alpha[t][j] += alpha[t-1][i] * a[i][j];
cannam@19 440 alpha[t][j] *= b[t][j];
cannam@19 441
cannam@19 442 c[t] += alpha[t][j];
cannam@19 443 }
cannam@19 444
cannam@19 445 /*
cannam@19 446 if (c[t] == 0)
cannam@19 447 {
cannam@19 448 printf("c[%d] = 0, going to blow up so exiting\n", t);
cannam@19 449 for (i = 0; i < N; i++)
cannam@19 450 if (b[t][i] == 0)
cannam@19 451 fprintf(stderr, "b[%d][%d] was zero\n", t, i);
cannam@19 452 fprintf(stderr, "x[t] was \n");
cannam@19 453 for (i = 0; i < L; i++)
cannam@19 454 fprintf(stderr, "%f ", x[t][i]);
cannam@19 455 fprintf(stderr, "\n\n");
cannam@19 456 exit(-1);
cannam@19 457 }
cannam@19 458 */
cannam@19 459
cannam@19 460 c[t] = 1 / c[t];
cannam@19 461 for (j = 0; j < N; j++)
cannam@19 462 alpha[t][j] *= c[t];
cannam@19 463
cannam@19 464 //printf("c[%d] = %e\n", t, c[t]);
cannam@19 465
cannam@19 466 *loglik -= log(c[t]);
cannam@19 467 }
cannam@19 468
cannam@19 469 /* calculate backwards probs using same coefficients */
cannam@19 470 for (i = 0; i < N; i++)
cannam@19 471 beta[T-1][i] = 1;
cannam@19 472 t = T - 1;
cannam@19 473 while (1)
cannam@19 474 {
cannam@19 475 for (i = 0; i < N; i++)
cannam@19 476 beta[t][i] *= c[t];
cannam@19 477
cannam@19 478 if (t == 0)
cannam@19 479 break;
cannam@19 480
cannam@19 481 for (i = 0; i < N; i++)
cannam@19 482 {
cannam@19 483 beta[t-1][i] = 0;
cannam@19 484 for (j = 0; j < N; j++)
cannam@19 485 beta[t-1][i] += a[i][j] * b[t][j] * beta[t][j];
cannam@19 486 }
cannam@19 487
cannam@19 488 t--;
cannam@19 489 }
cannam@19 490
cannam@19 491 /*
cannam@19 492 printf("alpha:\n");
cannam@19 493 for (t = 0; t < T; t++)
cannam@19 494 {
cannam@19 495 for (i = 0; i < N; i++)
cannam@19 496 printf("%4.4e\t\t", alpha[t][i]);
cannam@19 497 printf("\n");
cannam@19 498 }
cannam@19 499 printf("\n\n");printf("beta:\n");
cannam@19 500 for (t = 0; t < T; t++)
cannam@19 501 {
cannam@19 502 for (i = 0; i < N; i++)
cannam@19 503 printf("%4.4e\t\t", beta[t][i]);
cannam@19 504 printf("\n");
cannam@19 505 }
cannam@19 506 printf("\n\n");
cannam@19 507 */
cannam@19 508
cannam@19 509 /* calculate posterior probs */
cannam@19 510 double tot;
cannam@19 511 for (t = 0; t < T; t++)
cannam@19 512 {
cannam@19 513 tot = 0;
cannam@19 514 for (i = 0; i < N; i++)
cannam@19 515 {
cannam@19 516 gamma[t][i] = alpha[t][i] * beta[t][i];
cannam@19 517 tot += gamma[t][i];
cannam@19 518 }
cannam@19 519 for (i = 0; i < N; i++)
cannam@19 520 {
cannam@19 521 gamma[t][i] /= tot;
cannam@19 522
cannam@19 523 //printf("gamma[%d][%d] = %f\n", t, i, gamma[t][i]);
cannam@19 524 }
cannam@19 525 }
cannam@19 526
cannam@19 527 for (t = 0; t < T-1; t++)
cannam@19 528 {
cannam@19 529 tot = 0;
cannam@19 530 for (i = 0; i < N; i++)
cannam@19 531 {
cannam@19 532 for (j = 0; j < N; j++)
cannam@19 533 {
cannam@19 534 xi[t][i][j] = alpha[t][i] * a[i][j] * b[t+1][j] * beta[t+1][j];
cannam@19 535 tot += xi[t][i][j];
cannam@19 536 }
cannam@19 537 }
cannam@19 538 for (i = 0; i < N; i++)
cannam@19 539 for (j = 0; j < N; j++)
cannam@19 540 xi[t][i][j] /= tot;
cannam@19 541 }
cannam@19 542
cannam@19 543 /*
cannam@19 544 // CHECK - fine
cannam@19 545 // gamma[t][i] = \sum_j{xi[t][i][j]}
cannam@19 546 tot = 0;
cannam@19 547 for (j = 0; j < N; j++)
cannam@19 548 tot += xi[3][1][j];
cannam@19 549 printf("gamma[3][1] = %f, sum_j(xi[3][1][j]) = %f\n", gamma[3][1], tot);
cannam@19 550 */
cannam@19 551
cannam@19 552 for (t = 0; t < T; t++)
cannam@19 553 {
cannam@19 554 free(alpha[t]);
cannam@19 555 free(beta[t]);
cannam@19 556 }
cannam@19 557 free(alpha);
cannam@19 558 free(beta);
cannam@19 559 free(c);
cannam@19 560 }
cannam@19 561
cannam@19 562 void viterbi_decode(double** x, int T, model_t* model, int* q)
cannam@19 563 {
cannam@19 564 int i, j, t;
cannam@19 565 double max;
cannam@19 566
cannam@19 567 int N = model->N;
cannam@19 568 int L = model->L;
cannam@19 569 double* p0 = model->p0;
cannam@19 570 double** a = model->a;
cannam@19 571 double** mu = model->mu;
cannam@19 572 double** cov = model->cov;
cannam@19 573
cannam@19 574 /* inverse covariance and its determinant */
cannam@19 575 double** icov = (double**) malloc(L*sizeof(double*));
cannam@19 576 for (i = 0; i < L; i++)
cannam@19 577 icov[i] = (double*) malloc(L*sizeof(double));
cannam@19 578 double detcov;
cannam@19 579
cannam@19 580 double** logb = (double**) malloc(T*sizeof(double*));
cannam@19 581 double** phi = (double**) malloc(T*sizeof(double*));
cannam@19 582 int** psi = (int**) malloc(T*sizeof(int*));
cannam@19 583 for (t = 0; t < T; t++)
cannam@19 584 {
cannam@19 585 logb[t] = (double*) malloc(N*sizeof(double));
cannam@19 586 phi[t] = (double*) malloc(N*sizeof(double));
cannam@19 587 psi[t] = (int*) malloc(N*sizeof(int));
cannam@19 588 }
cannam@19 589
cannam@19 590 /* temporary memory */
cannam@19 591 double* gauss_y = (double*) malloc(L*sizeof(double));
cannam@19 592 double* gauss_z = (double*) malloc(L*sizeof(double));
cannam@19 593
cannam@19 594 /* calculate observation logprobs */
cannam@19 595 invert(cov, L, icov, &detcov);
cannam@19 596 for (t = 0; t < T; t++)
cannam@19 597 for (i = 0; i < N; i++)
cannam@19 598 logb[t][i] = loggauss(x[t], L, mu[i], icov, detcov, gauss_y, gauss_z);
cannam@19 599
cannam@19 600 /* initialise */
cannam@19 601 for (i = 0; i < N; i++)
cannam@19 602 {
cannam@19 603 phi[0][i] = log(p0[i]) + logb[0][i];
cannam@19 604 psi[0][i] = 0;
cannam@19 605 }
cannam@19 606
cannam@19 607 for (t = 1; t < T; t++)
cannam@19 608 {
cannam@19 609 for (j = 0; j < N; j++)
cannam@19 610 {
cannam@19 611 max = -1000000; // TODO: what should this be?? = smallest possible sumlogprob
cannam@19 612 psi[t][j] = 0;
cannam@19 613 for (i = 0; i < N; i++)
cannam@19 614 {
cannam@19 615 if (phi[t-1][i] + log(a[i][j]) > max)
cannam@19 616 {
cannam@19 617 max = phi[t-1][i] + log(a[i][j]);
cannam@19 618 phi[t][j] = max + logb[t][j];
cannam@19 619 psi[t][j] = i;
cannam@19 620 }
cannam@19 621 }
cannam@19 622 }
cannam@19 623 }
cannam@19 624
cannam@19 625 /* find maximising state at time T-1 */
cannam@19 626 max = phi[T-1][0];
cannam@19 627 q[T-1] = 0;
cannam@19 628 for (i = 1; i < N; i++)
cannam@19 629 {
cannam@19 630 if (phi[T-1][i] > max)
cannam@19 631 {
cannam@19 632 max = phi[T-1][i];
cannam@19 633 q[T-1] = i;
cannam@19 634 }
cannam@19 635 }
cannam@19 636
cannam@19 637
cannam@19 638 /* track back */
cannam@19 639 t = T - 2;
cannam@19 640 while (t >= 0)
cannam@19 641 {
cannam@19 642 q[t] = psi[t+1][q[t+1]];
cannam@19 643 t--;
cannam@19 644 }
cannam@19 645
cannam@19 646 /* de-allocate memory */
cannam@19 647 for (i = 0; i < L; i++)
cannam@19 648 free(icov[i]);
cannam@19 649 free(icov);
cannam@19 650 for (t = 0; t < T; t++)
cannam@19 651 {
cannam@19 652 free(logb[t]);
cannam@19 653 free(phi[t]);
cannam@19 654 free(psi[t]);
cannam@19 655 }
cannam@19 656 free(logb);
cannam@19 657 free(phi);
cannam@19 658 free(psi);
cannam@19 659
cannam@19 660 free(gauss_y);
cannam@19 661 free(gauss_z);
cannam@19 662 }
cannam@19 663
cannam@19 664 /* invert matrix and calculate determinant using LAPACK */
cannam@19 665 void invert(double** cov, int L, double** icov, double* detcov)
cannam@19 666 {
cannam@19 667 /* copy square matrix into a vector in column-major order */
cannam@19 668 double* a = (double*) malloc(L*L*sizeof(double));
cannam@19 669 int i, j;
cannam@19 670 for(j=0; j < L; j++)
cannam@19 671 for (i=0; i < L; i++)
cannam@19 672 a[j*L+i] = cov[i][j];
cannam@19 673
cannam@19 674 long M = (long) L;
cannam@19 675 long* ipiv = (long*) malloc(L*L*sizeof(int));
cannam@19 676 long ret;
cannam@19 677
cannam@19 678 /* LU decomposition */
cannam@19 679 ret = dgetrf_(&M, &M, a, &M, ipiv, &ret); /* ret should be zero, negative if cov is singular */
cannam@19 680 if (ret < 0)
cannam@19 681 {
cannam@19 682 fprintf(stderr, "Covariance matrix was singular, couldn't invert\n");
cannam@19 683 exit(-1);
cannam@19 684 }
cannam@19 685
cannam@19 686 /* find determinant */
cannam@19 687 double det = 1;
cannam@19 688 for(i = 0; i < L; i++)
cannam@19 689 det *= a[i*L+i];
cannam@19 690 // TODO: get this to work!!! If detcov < 0 then cov is bad anyway...
cannam@19 691 /*
cannam@19 692 int sign = 1;
cannam@19 693 for (i = 0; i < L; i++)
cannam@19 694 if (ipiv[i] != i)
cannam@19 695 sign = -sign;
cannam@19 696 det *= sign;
cannam@19 697 */
cannam@19 698 if (det < 0)
cannam@19 699 det = -det;
cannam@19 700 *detcov = det;
cannam@19 701
cannam@19 702 /* allocate required working storage */
cannam@19 703 long lwork = -1;
cannam@19 704 double lwbest;
cannam@19 705 dgetri_(&M, a, &M, ipiv, &lwbest, &lwork, &ret);
cannam@19 706 lwork = (long) lwbest;
cannam@19 707 double* work = (double*) malloc(lwork*sizeof(double));
cannam@19 708
cannam@19 709 /* find inverse */
cannam@19 710 dgetri_(&M, a, &M, ipiv, work, &lwork, &ret);
cannam@19 711
cannam@19 712 for(j=0; j < L; j++)
cannam@19 713 for (i=0; i < L; i++)
cannam@19 714 icov[i][j] = a[j*L+i];
cannam@19 715
cannam@19 716 free(work);
cannam@19 717 free(a);
cannam@19 718 }
cannam@19 719
cannam@19 720 /* probability of multivariate Gaussian given mean, inverse and determinant of covariance */
cannam@19 721 double gauss(double* x, int L, double* mu, double** icov, double detcov, double* y, double* z)
cannam@19 722 {
cannam@19 723 int i, j;
cannam@19 724 double s = 0;
cannam@19 725 for (i = 0; i < L; i++)
cannam@19 726 y[i] = x[i] - mu[i];
cannam@19 727 for (i = 0; i < L; i++)
cannam@19 728 {
cannam@19 729 //z[i] = 0;
cannam@19 730 //for (j = 0; j < L; j++)
cannam@19 731 // z[i] += icov[i][j] * y[j];
cannam@19 732 z[i] = cblas_ddot(L, &icov[i][0], 1, y, 1);
cannam@19 733 }
cannam@19 734 s = cblas_ddot(L, z, 1, y, 1);
cannam@19 735 //for (i = 0; i < L; i++)
cannam@19 736 // s += z[i] * y[i];
cannam@19 737
cannam@19 738 return exp(-s/2.0) / (pow(2*PI, L/2.0) * sqrt(detcov));
cannam@19 739 }
cannam@19 740
cannam@19 741 /* log probability of multivariate Gaussian given mean, inverse and determinant of covariance */
cannam@19 742 double loggauss(double* x, int L, double* mu, double** icov, double detcov, double* y, double* z)
cannam@19 743 {
cannam@19 744 int i, j;
cannam@19 745 double s = 0;
cannam@19 746 double ret;
cannam@19 747 for (i = 0; i < L; i++)
cannam@19 748 y[i] = x[i] - mu[i];
cannam@19 749 for (i = 0; i < L; i++)
cannam@19 750 {
cannam@19 751 //z[i] = 0;
cannam@19 752 //for (j = 0; j < L; j++)
cannam@19 753 // z[i] += icov[i][j] * y[j];
cannam@19 754 z[i] = cblas_ddot(L, &icov[i][0], 1, y, 1);
cannam@19 755 }
cannam@19 756 s = cblas_ddot(L, z, 1, y, 1);
cannam@19 757 //for (i = 0; i < L; i++)
cannam@19 758 // s += z[i] * y[i];
cannam@19 759
cannam@19 760 ret = -0.5 * (s + L * log(2*PI) + log(detcov));
cannam@19 761
cannam@19 762 /*
cannam@19 763 // TEST
cannam@19 764 if (isinf(ret) > 0)
cannam@19 765 printf("loggauss returning infinity\n");
cannam@19 766 if (isinf(ret) < 0)
cannam@19 767 printf("loggauss returning -infinity\n");
cannam@19 768 if (isnan(ret))
cannam@19 769 printf("loggauss returning nan\n");
cannam@19 770 */
cannam@19 771
cannam@19 772 return ret;
cannam@19 773 }
cannam@19 774
cannam@19 775 void hmm_print(model_t* model)
cannam@19 776 {
cannam@19 777 int i, j;
cannam@19 778 printf("p0:\n");
cannam@19 779 for (i = 0; i < model->N; i++)
cannam@19 780 printf("%f ", model->p0[i]);
cannam@19 781 printf("\n\n");
cannam@19 782 printf("a:\n");
cannam@19 783 for (i = 0; i < model->N; i++)
cannam@19 784 {
cannam@19 785 for (j = 0; j < model->N; j++)
cannam@19 786 printf("%f ", model->a[i][j]);
cannam@19 787 printf("\n");
cannam@19 788 }
cannam@19 789 printf("\n\n");
cannam@19 790 printf("mu:\n");
cannam@19 791 for (i = 0; i < model->N; i++)
cannam@19 792 {
cannam@19 793 for (j = 0; j < model->L; j++)
cannam@19 794 printf("%f ", model->mu[i][j]);
cannam@19 795 printf("\n");
cannam@19 796 }
cannam@19 797 printf("\n\n");
cannam@19 798 printf("cov:\n");
cannam@19 799 for (i = 0; i < model->L; i++)
cannam@19 800 {
cannam@19 801 for (j = 0; j < model->L; j++)
cannam@19 802 printf("%f ", model->cov[i][j]);
cannam@19 803 printf("\n");
cannam@19 804 }
cannam@19 805 printf("\n\n");
cannam@19 806 }
cannam@19 807
cannam@19 808