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