annotate hmm/HMM.c @ 242:6060110dc3c6

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