c@244: /* c@244: * hmm.c c@244: * soundbite c@244: * c@244: * Created by Mark Levy on 12/02/2006. c@244: * Copyright 2006 Centre for Digital Music, Queen Mary, University of London. All rights reserved. c@244: * c@244: */ c@244: c@244: #include c@244: #include c@244: #include c@244: #include c@244: #include /* to seed random number generator */ c@269: c@244: #include /* LAPACK for matrix inversion */ c@269: c@269: #ifdef ATLAS_ORDER c@269: #define HAVE_ATLAS 1 c@269: #endif c@269: c@269: #ifdef HAVE_ATLAS c@269: // Using ATLAS C interface to LAPACK c@269: #define dgetrf_(m, n, a, lda, ipiv, info) \ c@269: clapack_dgetrf(CblasColMajor, *m, *n, a, *lda, ipiv) c@269: #define dgetri_(n, a, lda, ipiv, work, lwork, info) \ c@269: clapack_dgetri(CblasColMajor, *n, a, *lda, ipiv) c@269: #endif c@269: c@244: #ifdef _MAC_OS_X c@244: #include c@244: #else c@244: #include /* BLAS for matrix multiplication */ c@244: #endif c@244: c@244: #include "hmm.h" c@244: c@244: model_t* hmm_init(double** x, int T, int L, int N) c@244: { c@244: int i, j, d, e, t; c@244: double s, ss; c@244: c@244: model_t* model; c@244: model = (model_t*) malloc(sizeof(model_t)); c@244: model->N = N; c@244: model->L = L; c@244: model->p0 = (double*) malloc(N*sizeof(double)); c@244: model->a = (double**) malloc(N*sizeof(double*)); c@244: model->mu = (double**) malloc(N*sizeof(double*)); c@244: for (i = 0; i < N; i++) c@244: { c@244: model->a[i] = (double*) malloc(N*sizeof(double)); c@244: model->mu[i] = (double*) malloc(L*sizeof(double)); c@244: } c@244: model->cov = (double**) malloc(L*sizeof(double*)); c@244: for (i = 0; i < L; i++) c@244: model->cov[i] = (double*) malloc(L*sizeof(double)); c@244: c@244: srand(time(0)); c@244: double* global_mean = (double*) malloc(L*sizeof(double)); c@244: c@244: /* find global mean */ c@244: for (d = 0; d < L; d++) c@244: { c@244: global_mean[d] = 0; c@244: for (t = 0; t < T; t++) c@244: global_mean[d] += x[t][d]; c@244: global_mean[d] /= T; c@244: } c@244: c@244: /* calculate global diagonal covariance */ c@244: for (d = 0; d < L; d++) c@244: { c@244: for (e = 0; e < L; e++) c@244: model->cov[d][e] = 0; c@244: for (t = 0; t < T; t++) c@244: model->cov[d][d] += (x[t][d] - global_mean[d]) * (x[t][d] - global_mean[d]); c@244: model->cov[d][d] /= T-1; c@244: } c@244: c@244: /* set all means close to global mean */ c@244: for (i = 0; i < N; i++) c@244: { c@244: for (d = 0; d < L; d++) c@244: { c@244: /* add some random noise related to covariance */ c@244: /* ideally the random number would be Gaussian(0,1), as a hack we make it uniform on [-0.25,0.25] */ c@244: model->mu[i][d] = global_mean[d] + (0.5 * rand() / (double) RAND_MAX - 0.25) * sqrt(model->cov[d][d]); c@244: } c@244: } c@244: c@244: /* random intial and transition probs */ c@244: s = 0; c@244: for (i = 0; i < N; i++) c@244: { c@244: model->p0[i] = 1 + rand() / (double) RAND_MAX; c@244: s += model->p0[i]; c@244: ss = 0; c@244: for (j = 0; j < N; j++) c@244: { c@244: model->a[i][j] = 1 + rand() / (double) RAND_MAX; c@244: ss += model->a[i][j]; c@244: } c@244: for (j = 0; j < N; j++) c@244: { c@244: model->a[i][j] /= ss; c@244: } c@244: } c@244: for (i = 0; i < N; i++) c@244: model->p0[i] /= s; c@244: c@244: free(global_mean); c@244: c@244: return model; c@244: } c@244: c@244: void hmm_close(model_t* model) c@244: { c@244: int i; c@244: c@244: for (i = 0; i < model->N; i++) c@244: { c@244: free(model->a[i]); c@244: free(model->mu[i]); c@244: } c@244: free(model->a); c@244: free(model->mu); c@244: for (i = 0; i < model->L; i++) c@244: free(model->cov[i]); c@244: free(model->cov); c@244: free(model); c@244: } c@244: c@244: void hmm_train(double** x, int T, model_t* model) c@244: { c@244: int i, t; c@244: double loglik; /* overall log-likelihood at each iteration */ c@244: c@244: int N = model->N; c@244: int L = model->L; c@244: double* p0 = model->p0; c@244: double** a = model->a; c@244: double** mu = model->mu; c@244: double** cov = model->cov; c@244: c@244: /* allocate memory */ c@244: double** gamma = (double**) malloc(T*sizeof(double*)); c@244: double*** xi = (double***) malloc(T*sizeof(double**)); c@244: for (t = 0; t < T; t++) c@244: { c@244: gamma[t] = (double*) malloc(N*sizeof(double)); c@244: xi[t] = (double**) malloc(N*sizeof(double*)); c@244: for (i = 0; i < N; i++) c@244: xi[t][i] = (double*) malloc(N*sizeof(double)); c@244: } c@244: c@244: /* temporary memory */ c@244: double* gauss_y = (double*) malloc(L*sizeof(double)); c@244: double* gauss_z = (double*) malloc(L*sizeof(double)); c@244: c@244: /* obs probs P(j|{x}) */ c@244: double** b = (double**) malloc(T*sizeof(double*)); c@244: for (t = 0; t < T; t++) c@244: b[t] = (double*) malloc(N*sizeof(double)); c@244: c@244: /* inverse covariance and its determinant */ c@244: double** icov = (double**) malloc(L*sizeof(double*)); c@244: for (i = 0; i < L; i++) c@244: icov[i] = (double*) malloc(L*sizeof(double)); c@244: double detcov; c@244: c@244: double thresh = 0.0001; c@244: int niter = 50; c@244: int iter = 0; c@244: double loglik1, loglik2; c@255: int foundnan = 0; c@255: c@255: while (iter < niter && !foundnan && !(iter > 1 && (loglik - loglik1) < thresh * (loglik1 - loglik2))) c@244: { c@244: ++iter; c@244: c@244: fprintf(stderr, "calculating obsprobs...\n"); c@244: fflush(stderr); c@244: c@244: /* precalculate obs probs */ c@244: invert(cov, L, icov, &detcov); c@244: c@244: for (t = 0; t < T; t++) c@244: { c@244: //int allzero = 1; c@244: for (i = 0; i < N; i++) c@244: { c@244: b[t][i] = exp(loggauss(x[t], L, mu[i], icov, detcov, gauss_y, gauss_z)); c@244: c@244: //if (b[t][i] != 0) c@244: // allzero = 0; c@244: } c@244: /* c@244: if (allzero) c@244: { c@244: printf("all the b[t][i] were zero for t = %d, correcting...\n", t); c@244: for (i = 0; i < N; i++) c@244: { c@244: b[t][i] = 0.00001; c@244: } c@244: } c@244: */ c@244: } c@244: c@244: fprintf(stderr, "forwards-backwards...\n"); c@244: fflush(stderr); c@244: c@244: forward_backwards(xi, gamma, &loglik, &loglik1, &loglik2, iter, N, T, p0, a, b); c@244: c@244: fprintf(stderr, "iteration %d: loglik = %f\n", iter, loglik); c@244: fprintf(stderr, "re-estimation...\n"); c@244: fflush(stderr); c@255: c@255: if (isnan(loglik)) { c@255: foundnan = 1; c@255: continue; c@255: } c@244: c@244: baum_welch(p0, a, mu, cov, N, T, L, x, xi, gamma); c@244: c@244: /* c@244: printf("a:\n"); c@244: for (i = 0; i < model->N; i++) c@244: { c@244: for (j = 0; j < model->N; j++) c@244: printf("%f ", model->a[i][j]); c@244: printf("\n"); c@244: } c@244: printf("\n\n"); c@244: */ c@244: //hmm_print(model); c@244: } c@244: c@244: /* deallocate memory */ c@244: for (t = 0; t < T; t++) c@244: { c@244: free(gamma[t]); c@244: free(b[t]); c@244: for (i = 0; i < N; i++) c@244: free(xi[t][i]); c@244: free(xi[t]); c@244: } c@244: free(gamma); c@244: free(xi); c@244: free(b); c@244: c@244: for (i = 0; i < L; i++) c@244: free(icov[i]); c@244: free(icov); c@244: c@244: free(gauss_y); c@244: free(gauss_z); c@244: } c@244: c@244: void mlss_reestimate(double* p0, double** a, double** mu, double** cov, int N, int T, int L, int* q, double** x) c@244: { c@244: /* fit a single Gaussian to observations in each state */ c@244: c@244: /* calculate the mean observation in each state */ c@244: c@244: /* calculate the overall covariance */ c@244: c@244: /* count transitions */ c@244: c@244: /* estimate initial probs from transitions (???) */ c@244: } c@244: c@244: 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: { c@244: int i, j, t; c@244: c@244: double* sum_gamma = (double*) malloc(N*sizeof(double)); c@244: c@244: /* temporary memory */ c@244: double* u = (double*) malloc(L*L*sizeof(double)); c@244: double* yy = (double*) malloc(T*L*sizeof(double)); c@244: double* yy2 = (double*) malloc(T*L*sizeof(double)); c@244: c@244: /* re-estimate transition probs */ c@244: for (i = 0; i < N; i++) c@244: { c@244: sum_gamma[i] = 0; c@244: for (t = 0; t < T-1; t++) c@244: sum_gamma[i] += gamma[t][i]; c@244: } c@244: c@244: for (i = 0; i < N; i++) c@244: { c@244: if (sum_gamma[i] == 0) c@244: { c@244: fprintf(stderr, "sum_gamma[%d] was zero...\n", i); c@244: } c@244: //double s = 0; c@244: for (j = 0; j < N; j++) c@244: { c@244: a[i][j] = 0; c@255: if (sum_gamma[i] == 0.) continue; c@244: for (t = 0; t < T-1; t++) c@244: a[i][j] += xi[t][i][j]; c@244: //s += a[i][j]; c@244: a[i][j] /= sum_gamma[i]; c@244: } c@244: /* c@244: for (j = 0; j < N; j++) c@244: { c@244: a[i][j] /= s; c@244: } c@244: */ c@244: } c@244: c@244: /* NB: now we need to sum gamma over all t */ c@244: for (i = 0; i < N; i++) c@244: sum_gamma[i] += gamma[T-1][i]; c@244: c@244: /* re-estimate initial probs */ c@244: for (i = 0; i < N; i++) c@244: p0[i] = gamma[0][i]; c@244: c@244: /* re-estimate covariance */ c@244: int d, e; c@244: double sum_sum_gamma = 0; c@244: for (i = 0; i < N; i++) c@244: sum_sum_gamma += sum_gamma[i]; c@244: c@244: /* c@244: for (d = 0; d < L; d++) c@244: { c@244: for (e = d; e < L; e++) c@244: { c@244: cov[d][e] = 0; c@244: for (t = 0; t < T; t++) c@244: for (j = 0; j < N; j++) c@244: cov[d][e] += gamma[t][j] * (x[t][d] - mu[j][d]) * (x[t][e] - mu[j][e]); c@244: c@244: cov[d][e] /= sum_sum_gamma; c@244: c@244: if (isnan(cov[d][e])) c@244: { c@244: printf("cov[%d][%d] was nan\n", d, e); c@244: for (j = 0; j < N; j++) c@244: for (i = 0; i < L; i++) c@244: if (isnan(mu[j][i])) c@244: printf("mu[%d][%d] was nan\n", j, i); c@244: for (t = 0; t < T; t++) c@244: for (j = 0; j < N; j++) c@244: if (isnan(gamma[t][j])) c@244: printf("gamma[%d][%d] was nan\n", t, j); c@244: exit(-1); c@244: } c@244: } c@244: } c@244: for (d = 0; d < L; d++) c@244: for (e = 0; e < d; e++) c@244: cov[d][e] = cov[e][d]; c@244: */ c@244: c@244: /* using BLAS */ c@244: for (d = 0; d < L; d++) c@244: for (e = 0; e < L; e++) c@244: cov[d][e] = 0; c@244: c@244: for (j = 0; j < N; j++) c@244: { c@244: for (d = 0; d < L; d++) c@244: for (t = 0; t < T; t++) c@244: { c@244: yy[d*T+t] = x[t][d] - mu[j][d]; c@244: yy2[d*T+t] = gamma[t][j] * (x[t][d] - mu[j][d]); c@244: } c@244: c@244: cblas_dgemm(CblasColMajor, CblasTrans, CblasNoTrans, L, L, T, 1.0, yy, T, yy2, T, 0, u, L); c@244: c@244: for (e = 0; e < L; e++) c@244: for (d = 0; d < L; d++) c@244: cov[d][e] += u[e*L+d]; c@244: } c@244: c@244: for (d = 0; d < L; d++) c@244: for (e = 0; e < L; e++) c@244: cov[d][e] /= T; /* sum_sum_gamma; */ c@244: c@244: //printf("sum_sum_gamma = %f\n", sum_sum_gamma); /* fine, = T IS THIS ALWAYS TRUE with pooled cov?? */ c@244: c@244: /* re-estimate means */ c@244: for (j = 0; j < N; j++) c@244: { c@244: for (d = 0; d < L; d++) c@244: { c@244: mu[j][d] = 0; c@244: for (t = 0; t < T; t++) c@244: mu[j][d] += gamma[t][j] * x[t][d]; c@244: mu[j][d] /= sum_gamma[j]; c@244: } c@244: } c@244: c@244: /* deallocate memory */ c@244: free(sum_gamma); c@244: free(yy); c@244: free(yy2); c@244: free(u); c@244: } c@244: c@244: 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: { c@244: /* forwards-backwards with scaling */ c@244: int i, j, t; c@244: c@244: double** alpha = (double**) malloc(T*sizeof(double*)); c@244: double** beta = (double**) malloc(T*sizeof(double*)); c@244: for (t = 0; t < T; t++) c@244: { c@244: alpha[t] = (double*) malloc(N*sizeof(double)); c@244: beta[t] = (double*) malloc(N*sizeof(double)); c@244: } c@244: c@244: /* scaling coefficients */ c@244: double* c = (double*) malloc(T*sizeof(double)); c@244: c@244: /* calculate forward probs and scale coefficients */ c@244: c[0] = 0; c@244: for (i = 0; i < N; i++) c@244: { c@244: alpha[0][i] = p0[i] * b[0][i]; c@244: c[0] += alpha[0][i]; c@244: c@244: //printf("p0[%d] = %f, b[0][%d] = %f\n", i, p0[i], i, b[0][i]); c@244: } c@244: c[0] = 1 / c[0]; c@244: for (i = 0; i < N; i++) c@244: { c@244: alpha[0][i] *= c[0]; c@244: c@244: //printf("alpha[0][%d] = %f\n", i, alpha[0][i]); /* OK agrees with Matlab */ c@244: } c@244: c@244: *loglik1 = *loglik; c@244: *loglik = -log(c[0]); c@244: if (iter == 2) c@244: *loglik2 = *loglik; c@244: c@244: for (t = 1; t < T; t++) c@244: { c@244: c[t] = 0; c@244: for (j = 0; j < N; j++) c@244: { c@244: alpha[t][j] = 0; c@244: for (i = 0; i < N; i++) c@244: alpha[t][j] += alpha[t-1][i] * a[i][j]; c@244: alpha[t][j] *= b[t][j]; c@244: c@244: c[t] += alpha[t][j]; c@244: } c@244: c@244: /* c@244: if (c[t] == 0) c@244: { c@244: printf("c[%d] = 0, going to blow up so exiting\n", t); c@244: for (i = 0; i < N; i++) c@244: if (b[t][i] == 0) c@244: fprintf(stderr, "b[%d][%d] was zero\n", t, i); c@244: fprintf(stderr, "x[t] was \n"); c@244: for (i = 0; i < L; i++) c@244: fprintf(stderr, "%f ", x[t][i]); c@244: fprintf(stderr, "\n\n"); c@244: exit(-1); c@244: } c@244: */ c@244: c@244: c[t] = 1 / c[t]; c@244: for (j = 0; j < N; j++) c@244: alpha[t][j] *= c[t]; c@244: c@244: //printf("c[%d] = %e\n", t, c[t]); c@244: c@244: *loglik -= log(c[t]); c@244: } c@244: c@244: /* calculate backwards probs using same coefficients */ c@244: for (i = 0; i < N; i++) c@244: beta[T-1][i] = 1; c@244: t = T - 1; c@244: while (1) c@244: { c@244: for (i = 0; i < N; i++) c@244: beta[t][i] *= c[t]; c@244: c@244: if (t == 0) c@244: break; c@244: c@244: for (i = 0; i < N; i++) c@244: { c@244: beta[t-1][i] = 0; c@244: for (j = 0; j < N; j++) c@244: beta[t-1][i] += a[i][j] * b[t][j] * beta[t][j]; c@244: } c@244: c@244: t--; c@244: } c@244: c@244: /* c@244: printf("alpha:\n"); c@244: for (t = 0; t < T; t++) c@244: { c@244: for (i = 0; i < N; i++) c@244: printf("%4.4e\t\t", alpha[t][i]); c@244: printf("\n"); c@244: } c@244: printf("\n\n");printf("beta:\n"); c@244: for (t = 0; t < T; t++) c@244: { c@244: for (i = 0; i < N; i++) c@244: printf("%4.4e\t\t", beta[t][i]); c@244: printf("\n"); c@244: } c@244: printf("\n\n"); c@244: */ c@244: c@244: /* calculate posterior probs */ c@244: double tot; c@244: for (t = 0; t < T; t++) c@244: { c@244: tot = 0; c@244: for (i = 0; i < N; i++) c@244: { c@244: gamma[t][i] = alpha[t][i] * beta[t][i]; c@244: tot += gamma[t][i]; c@244: } c@244: for (i = 0; i < N; i++) c@244: { c@244: gamma[t][i] /= tot; c@244: c@244: //printf("gamma[%d][%d] = %f\n", t, i, gamma[t][i]); c@244: } c@244: } c@244: c@244: for (t = 0; t < T-1; t++) c@244: { c@244: tot = 0; c@244: for (i = 0; i < N; i++) c@244: { c@244: for (j = 0; j < N; j++) c@244: { c@244: xi[t][i][j] = alpha[t][i] * a[i][j] * b[t+1][j] * beta[t+1][j]; c@244: tot += xi[t][i][j]; c@244: } c@244: } c@244: for (i = 0; i < N; i++) c@244: for (j = 0; j < N; j++) c@244: xi[t][i][j] /= tot; c@244: } c@244: c@244: /* c@244: // CHECK - fine c@244: // gamma[t][i] = \sum_j{xi[t][i][j]} c@244: tot = 0; c@244: for (j = 0; j < N; j++) c@244: tot += xi[3][1][j]; c@244: printf("gamma[3][1] = %f, sum_j(xi[3][1][j]) = %f\n", gamma[3][1], tot); c@244: */ c@244: c@244: for (t = 0; t < T; t++) c@244: { c@244: free(alpha[t]); c@244: free(beta[t]); c@244: } c@244: free(alpha); c@244: free(beta); c@244: free(c); c@244: } c@244: c@244: void viterbi_decode(double** x, int T, model_t* model, int* q) c@244: { c@244: int i, j, t; c@244: double max; c@273: int havemax; c@244: c@244: int N = model->N; c@244: int L = model->L; c@244: double* p0 = model->p0; c@244: double** a = model->a; c@244: double** mu = model->mu; c@244: double** cov = model->cov; c@244: c@244: /* inverse covariance and its determinant */ c@244: double** icov = (double**) malloc(L*sizeof(double*)); c@244: for (i = 0; i < L; i++) c@244: icov[i] = (double*) malloc(L*sizeof(double)); c@244: double detcov; c@244: c@244: double** logb = (double**) malloc(T*sizeof(double*)); c@244: double** phi = (double**) malloc(T*sizeof(double*)); c@244: int** psi = (int**) malloc(T*sizeof(int*)); c@244: for (t = 0; t < T; t++) c@244: { c@244: logb[t] = (double*) malloc(N*sizeof(double)); c@244: phi[t] = (double*) malloc(N*sizeof(double)); c@244: psi[t] = (int*) malloc(N*sizeof(int)); c@244: } c@244: c@244: /* temporary memory */ c@244: double* gauss_y = (double*) malloc(L*sizeof(double)); c@244: double* gauss_z = (double*) malloc(L*sizeof(double)); c@244: c@244: /* calculate observation logprobs */ c@244: invert(cov, L, icov, &detcov); c@244: for (t = 0; t < T; t++) c@244: for (i = 0; i < N; i++) c@244: logb[t][i] = loggauss(x[t], L, mu[i], icov, detcov, gauss_y, gauss_z); c@244: c@244: /* initialise */ c@244: for (i = 0; i < N; i++) c@244: { c@244: phi[0][i] = log(p0[i]) + logb[0][i]; c@244: psi[0][i] = 0; c@244: } c@244: c@244: for (t = 1; t < T; t++) c@244: { c@244: for (j = 0; j < N; j++) c@244: { c@273: max = -1000000; c@273: havemax = 0; c@273: c@244: psi[t][j] = 0; c@244: for (i = 0; i < N; i++) c@244: { c@273: if (phi[t-1][i] + log(a[i][j]) > max || !havemax) c@244: { c@244: max = phi[t-1][i] + log(a[i][j]); c@244: phi[t][j] = max + logb[t][j]; c@244: psi[t][j] = i; c@273: havemax = 1; c@244: } c@244: } c@244: } c@244: } c@244: c@244: /* find maximising state at time T-1 */ c@244: max = phi[T-1][0]; c@244: q[T-1] = 0; c@244: for (i = 1; i < N; i++) c@244: { c@244: if (phi[T-1][i] > max) c@244: { c@244: max = phi[T-1][i]; c@244: q[T-1] = i; c@244: } c@244: } c@244: c@244: c@244: /* track back */ c@244: t = T - 2; c@244: while (t >= 0) c@244: { c@244: q[t] = psi[t+1][q[t+1]]; c@244: t--; c@244: } c@244: c@244: /* de-allocate memory */ c@244: for (i = 0; i < L; i++) c@244: free(icov[i]); c@244: free(icov); c@244: for (t = 0; t < T; t++) c@244: { c@244: free(logb[t]); c@244: free(phi[t]); c@244: free(psi[t]); c@244: } c@244: free(logb); c@244: free(phi); c@244: free(psi); c@244: c@244: free(gauss_y); c@244: free(gauss_z); c@244: } c@244: c@244: /* invert matrix and calculate determinant using LAPACK */ c@244: void invert(double** cov, int L, double** icov, double* detcov) c@244: { c@244: /* copy square matrix into a vector in column-major order */ c@244: double* a = (double*) malloc(L*L*sizeof(double)); c@244: int i, j; c@244: for(j=0; j < L; j++) c@244: for (i=0; i < L; i++) c@244: a[j*L+i] = cov[i][j]; c@244: c@269: int M = (int) L; c@269: int* ipiv = (int *) malloc(L*L*sizeof(int)); c@269: int ret; c@244: c@244: /* LU decomposition */ c@244: ret = dgetrf_(&M, &M, a, &M, ipiv, &ret); /* ret should be zero, negative if cov is singular */ c@244: if (ret < 0) c@244: { c@244: fprintf(stderr, "Covariance matrix was singular, couldn't invert\n"); c@244: exit(-1); c@244: } c@244: c@244: /* find determinant */ c@244: double det = 1; c@244: for(i = 0; i < L; i++) c@244: det *= a[i*L+i]; c@244: // TODO: get this to work!!! If detcov < 0 then cov is bad anyway... c@244: /* c@244: int sign = 1; c@244: for (i = 0; i < L; i++) c@244: if (ipiv[i] != i) c@244: sign = -sign; c@244: det *= sign; c@244: */ c@244: if (det < 0) c@244: det = -det; c@244: *detcov = det; c@244: c@244: /* allocate required working storage */ c@269: #ifndef HAVE_ATLAS c@269: int lwork = -1; c@269: double lwbest = 0.0; c@244: dgetri_(&M, a, &M, ipiv, &lwbest, &lwork, &ret); c@269: lwork = (int) lwbest; c@244: double* work = (double*) malloc(lwork*sizeof(double)); c@269: #endif c@244: c@244: /* find inverse */ c@244: dgetri_(&M, a, &M, ipiv, work, &lwork, &ret); c@269: c@244: for(j=0; j < L; j++) c@244: for (i=0; i < L; i++) c@244: icov[i][j] = a[j*L+i]; c@244: c@269: #ifndef HAVE_ATLAS c@244: free(work); c@269: #endif c@244: free(a); c@244: } c@244: c@244: /* probability of multivariate Gaussian given mean, inverse and determinant of covariance */ c@244: double gauss(double* x, int L, double* mu, double** icov, double detcov, double* y, double* z) c@244: { c@244: int i, j; c@244: double s = 0; c@244: for (i = 0; i < L; i++) c@244: y[i] = x[i] - mu[i]; c@244: for (i = 0; i < L; i++) c@244: { c@244: //z[i] = 0; c@244: //for (j = 0; j < L; j++) c@244: // z[i] += icov[i][j] * y[j]; c@244: z[i] = cblas_ddot(L, &icov[i][0], 1, y, 1); c@244: } c@244: s = cblas_ddot(L, z, 1, y, 1); c@244: //for (i = 0; i < L; i++) c@244: // s += z[i] * y[i]; c@244: c@244: return exp(-s/2.0) / (pow(2*PI, L/2.0) * sqrt(detcov)); c@244: } c@244: c@244: /* log probability of multivariate Gaussian given mean, inverse and determinant of covariance */ c@244: double loggauss(double* x, int L, double* mu, double** icov, double detcov, double* y, double* z) c@244: { c@244: int i, j; c@244: double s = 0; c@244: double ret; c@244: for (i = 0; i < L; i++) c@244: y[i] = x[i] - mu[i]; c@244: for (i = 0; i < L; i++) c@244: { c@244: //z[i] = 0; c@244: //for (j = 0; j < L; j++) c@244: // z[i] += icov[i][j] * y[j]; c@244: z[i] = cblas_ddot(L, &icov[i][0], 1, y, 1); c@244: } c@244: s = cblas_ddot(L, z, 1, y, 1); c@244: //for (i = 0; i < L; i++) c@244: // s += z[i] * y[i]; c@244: c@244: ret = -0.5 * (s + L * log(2*PI) + log(detcov)); c@244: c@244: /* c@244: // TEST c@244: if (isinf(ret) > 0) c@244: printf("loggauss returning infinity\n"); c@244: if (isinf(ret) < 0) c@244: printf("loggauss returning -infinity\n"); c@244: if (isnan(ret)) c@244: printf("loggauss returning nan\n"); c@244: */ c@244: c@244: return ret; c@244: } c@244: c@244: void hmm_print(model_t* model) c@244: { c@244: int i, j; c@244: printf("p0:\n"); c@244: for (i = 0; i < model->N; i++) c@244: printf("%f ", model->p0[i]); c@244: printf("\n\n"); c@244: printf("a:\n"); c@244: for (i = 0; i < model->N; i++) c@244: { c@244: for (j = 0; j < model->N; j++) c@244: printf("%f ", model->a[i][j]); c@244: printf("\n"); c@244: } c@244: printf("\n\n"); c@244: printf("mu:\n"); c@244: for (i = 0; i < model->N; i++) c@244: { c@244: for (j = 0; j < model->L; j++) c@244: printf("%f ", model->mu[i][j]); c@244: printf("\n"); c@244: } c@244: printf("\n\n"); c@244: printf("cov:\n"); c@244: for (i = 0; i < model->L; i++) c@244: { c@244: for (j = 0; j < model->L; j++) c@244: printf("%f ", model->cov[i][j]); c@244: printf("\n"); c@244: } c@244: printf("\n\n"); c@244: } c@244: c@244: