annotate hmm/hmm.c @ 189:e4a57215ddee

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