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