annotate hmm/hmm.c @ 269:a63c7b6191b5

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