annotate hmm/hmm.c @ 304:702ff8c08137

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