annotate hmm/hmm.c @ 308:25af9a1e4ec3

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