annotate hmm/hmm.c @ 483:fdaa63607c15

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