annotate hmm/HMM.c @ 18:8e90a56b4b5f

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