forked from OHBA-analysis/HMM-MAR
-
Notifications
You must be signed in to change notification settings - Fork 0
/
hmmdecode.m
382 lines (323 loc) · 11.7 KB
/
hmmdecode.m
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
function [Path,Xi] = hmmdecode(data,T,hmm,type,residuals,preproc)
%
% State time course and Viterbi decoding for hmm
% The algorithm is run for the whole data set, including those whose class
% was fixed. This means that the assignment for those can be different.
%
% INPUT
% data observations, either a struct with X (time series) and C (classes, optional)
% or just a matrix containing the time series
% T length of series
% hmm hmm data structure
% type 0, state time courses (default); 1, viterbi path
% residuals in case we train on residuals, the value of those (optional)
% preproc whether we should perform the preprocessing options with
% which the hmm model was trained; 1 by default.
%
% OUTPUT
% vpath (T x 1) maximum likelihood state sequence (type=1) OR
% vpath (T x K) state time courses
% Xi joint probability of past and future states conditioned on data
% (empty if Viterbi path is computed)
%
% Author: Diego Vidaurre, OHBA, University of Oxford
% to fix potential compatibility issues with previous versions
hmm = versCompatibilityFix(hmm);
if nargin<4 || isempty(type), type = 0; end
if nargin<5, residuals = []; end
if nargin<6 || isempty(preproc), preproc = 1; end
% if nargin<7 || isempty(grouping)
% if isfield(hmm.train,'grouping')
% grouping = hmm.train.grouping;
% else
% grouping = ones(length(T),1);
% end
% if size(grouping,1)==1, grouping = grouping'; end
% end
% if length(size(hmm.Dir_alpha))==3 && isempty(grouping)
% error('You must specify the grouping argument if the HMM was trained on different groups')
% elseif ~isempty(grouping)
% Q = length(unique(grouping));
% else
% Q = 1;
% end
stochastic_learn = isfield(hmm.train,'BIGNbatch') && hmm.train.BIGNbatch < length(T);
p = hmm.train.lowrank; do_HMM_pca = (p > 0);
if xor(iscell(data),iscell(T)), error('data and T must be cells, either both or none of them.'); end
if stochastic_learn
N = length(T);
if ~iscell(data)
dat = cell(N,1); TT = cell(N,1);
for i = 1:N
t = 1:T(i);
dat{i} = data(t,:); TT{i} = T(i);
try data(t,:) = [];
catch, error('The dimension of data does not correspond to T');
end
end
if ~isempty(data)
error('The dimension of data does not correspond to T');
end
data = dat; T = TT; clear dat TT
end
if nargin<2
Path = hmmsdecode(data,T,hmm,type);
else
[Path,Xi] = hmmsdecode(data,T,hmm,type);
end
return
else % data can be a cell or a matrix
if iscell(T)
for i = 1:length(T)
if size(T{i},1)==1, T{i} = T{i}'; end
end
if size(T,1)==1, T = T'; end
T = cell2mat(T);
end
checkdatacell;
N = length(T);
end
if preproc % Adjust the data if necessary
train = hmm.train;
checkdatacell;
data = data2struct(data,T,train);
% Standardise data and control for ackward trials
data = standardisedata(data,T,train.standardise);
% Filtering
if ~isempty(train.filter)
data = filterdata(data,T,train.Fs,train.filter);
end
% Detrend data
if train.detrend
data = detrenddata(data,T);
end
% Leakage correction
if train.leakagecorr ~= 0
data = leakcorr(data,T,train.leakagecorr);
end
% Hilbert envelope
if train.onpower
data = rawsignal2power(data,T);
end
% Leading Phase Eigenvectors
if train.leida
data = leadingPhEigenvector(data,T);
end
% pre-embedded PCA transform
if length(train.pca_spatial) > 1 || train.pca_spatial > 0
if isfield(train,'As')
data.X = bsxfun(@minus,data.X,mean(data.X));
data.X = data.X * train.As;
else
[train.As,data.X] = highdim_pca(data.X,T,train.pca_spatial);
end
end
% Embedding
if length(train.embeddedlags) > 1
[data,T] = embeddata(data,T,train.embeddedlags);
end
% PCA transform
if length(train.pca) > 1 || train.pca > 0
if isfield(train,'A')
data.X = bsxfun(@minus,data.X,mean(data.X));
data.X = data.X * train.A;
else
[train.A,data.X] = highdim_pca(data.X,T,train.pca,0,0,0,train.varimax);
end
% Standardise principal components and control for ackward trials
data = standardisedata(data,T,train.standardise_pc);
train.ndim = size(train.A,2);
train.S = ones(train.ndim);
orders = formorders(train.order,train.orderoffset,train.timelag,train.exptimelag);
train.Sind = formindexes(orders,train.S);
end
% Downsampling
if train.downsample > 0
[data,T] = downsampledata(data,T,train.downsample,train.Fs);
end
end
if type==0
[Path,~,Xi] = hsinference(data,T,hmm,residuals);
return
end
if isstruct(data)
if isfield(data,'C') && ~all(isnan(data.C(:)))
warning('Pre-specified state time courses will be ignored for Viterbi path calculation')
end
data = data.X;
end
Xi = [];
K = length(hmm.state);
if isempty(residuals) && ~do_HMM_pca
if ~isfield(hmm.train,'Sind')
orders = formorders(hmm.train.order,hmm.train.orderoffset,hmm.train.timelag,hmm.train.exptimelag);
hmm.train.Sind = formindexes(orders,hmm.train.S);
end
residuals = getresiduals(data,T,hmm.train.Sind,hmm.train.maxorder,hmm.train.order,...
hmm.train.orderoffset,hmm.train.timelag,hmm.train.exptimelag,hmm.train.zeromean);
end
if ~isfield(hmm,'P')
hmm = hmmhsinit(hmm);
end
order = hmm.train.maxorder;
if hmm.train.useParallel==1 && N>1
% to duplicate this code is really ugly but there doesn't seem to be
% any other way - more Matlab's fault than mine
Path = cell(N,1);
parfor n = 1:N
%if Q > 1
% i = grouping(n);
% P = hmm.P(:,:,i); Pi = hmm.Pi(:,i)';
%else
% P = hmm.P; Pi = hmm.Pi;
%end
% This causes error with the Parallel toolbox
P = hmm.P; Pi = hmm.Pi;
q_star = ones(T(n)-order,1);
scale=zeros(T(n),1);
alpha=zeros(T(n),K);
beta=zeros(T(n),K);
% Initialise Viterbi bits
delta=zeros(T(n),K);
psi=zeros(T(n),K);
if n==1, t0 = 0; s0 = 0;
else t0 = sum(T(1:n-1)); s0 = t0 - order*(n-1);
end
if do_HMM_pca
B = obslike(data(t0+1:t0+T(n),:),hmm,[]);
else
B = obslike(data(t0+1:t0+T(n),:),hmm,residuals(s0+1:s0+T(n)-order,:));
end
B(B<realmin) = realmin;
% Scaling for delta
dscale=zeros(T(n),1);
alpha(1+order,:)=Pi(:)'.*B(1+order,:);
scale(1+order)=sum(alpha(1+order,:));
alpha(1+order,:)=alpha(1+order,:)/(scale(1+order));
delta(1+order,:) = alpha(1+order,:); % Eq. 32(a) Rabiner (1989)
% Eq. 32(b) Psi already zero
for i=2+order:T(n)
alpha(i,:)=(alpha(i-1,:)*P).*B(i,:);
scale(i)=sum(alpha(i,:));
if scale(i)<realmin, scale(i) = realmin; end
alpha(i,:)=alpha(i,:)/(scale(i));
for k=1:K
v=delta(i-1,:).*P(:,k)';
mv=max(v);
delta(i,k)=mv*B(i,k); % Eq 33a Rabiner (1989)
fmv = find(v==mv);
if length(fmv) > 1
% no unique maximum - so pick one at random
tmp1=fmv;
tmp2=rand(length(tmp1),1);
[~,tmp4]=max(tmp2);
psi(i,k)=tmp4;
else
psi(i,k)=fmv; % ARGMAX; Eq 33b Rabiner (1989)
end
end
% SCALING FOR DELTA ????
dscale(i)=sum(delta(i,:));
if dscale(i)<realmin, dscale(i) = realmin; end
delta(i,:)=delta(i,:)/(dscale(i));
end
% Get beta values for single state decoding
beta(T(n),:)=ones(1,K)/scale(T(n));
for i=T(n)-1:-1:1+order
beta(i,:)=(beta(i+1,:).*B(i+1,:))*(P')/scale(i);
end
xi=zeros(T(n)-1-order,K*K);
for i=1+order:T(n)-1
t=P.*( alpha(i,:)' * (beta(i+1,:).*B(i+1,:)));
xi(i-order,:)=t(:)'/sum(t(:));
end
delta=delta(1+order:T(n),:);
psi=psi(1+order:T(n),:);
% Backtracking for Viterbi decoding
id = find(delta(T(n)-order,:)==max(delta(T(n)-order,:)));% Eq 34b Rabiner;
q_star(T(n)-order) = id(1);
for i=T(n)-1-order:-1:1
q_star(i) = psi(i+1,q_star(i+1));
end
Path{n} = single(q_star);
end
Path = cell2mat(Path);
else
Path = zeros(sum(T)-length(T)*order,1,'single');
tacc = 0;
for n=1:N
%if Q > 1
% i = grouping(n);
% P = hmm.P(:,:,i); Pi = hmm.Pi(:,i)';
%else
% P = hmm.P; Pi = hmm.Pi;
%end
P = hmm.P; Pi = hmm.Pi;
q_star = ones(T(n)-order,1);
alpha=zeros(T(n),K);
beta=zeros(T(n),K);
% Initialise Viterbi bits
delta=zeros(T(n),K);
psi=zeros(T(n),K);
if n==1, t0 = 0; s0 = 0;
else t0 = sum(T(1:n-1)); s0 = t0 - order*(n-1);
end
if do_HMM_pca
B = obslike(data(t0+1:t0+T(n),:),hmm,[]);
else
B = obslike(data(t0+1:t0+T(n),:),hmm,residuals(s0+1:s0+T(n)-order,:));
end
B(B<realmin) = realmin;
scale=zeros(T(n),1);
% Scaling for delta
dscale=zeros(T(n),1);
alpha(1+order,:)=Pi(:)'.*B(1+order,:);
scale(1+order)=sum(alpha(1+order,:));
alpha(1+order,:)=alpha(1+order,:)/(scale(1+order)+realmin);
delta(1+order,:) = alpha(1+order,:); % Eq. 32(a) Rabiner (1989)
% Eq. 32(b) Psi already zero
for i=2+order:T(n)
alpha(i,:)=(alpha(i-1,:)*P).*B(i,:);
scale(i)=sum(alpha(i,:));
alpha(i,:)=alpha(i,:)/(scale(i)+realmin);
for k=1:K
v=delta(i-1,:).*P(:,k)';
mv=max(v);
delta(i,k)=mv*B(i,k); % Eq 33a Rabiner (1989)
if length(find(v==mv)) > 1
% no unique maximum - so pick one at random
tmp1=find(v==mv);
tmp2=rand(length(tmp1),1);
[~,tmp4]=max(tmp2);
psi(i,k)=tmp4;
else
psi(i,k)=find(v==mv); % ARGMAX; Eq 33b Rabiner (1989)
end
end
% SCALING FOR DELTA ????
dscale(i)=sum(delta(i,:));
delta(i,:)=delta(i,:)/(dscale(i)+realmin);
end
% Get beta values for single state decoding
beta(T(n),:)=ones(1,K)/scale(T(n));
for i=T(n)-1:-1:1+order
beta(i,:)=(beta(i+1,:).*B(i+1,:))*(P')/scale(i);
end
xi=zeros(T(n)-1-order,K*K);
for i=1+order:T(n)-1
t=P.*( alpha(i,:)' * (beta(i+1,:).*B(i+1,:)));
xi(i-order,:)=t(:)'/sum(t(:));
end
delta=delta(1+order:T(n),:);
psi=psi(1+order:T(n),:);
% Backtracking for Viterbi decoding
id = find(delta(T(n)-order,:)==max(delta(T(n)-order,:)));% Eq 34b Rabiner;
q_star(T(n)-order) = id(1);
for i=T(n)-1-order:-1:1
q_star(i) = psi(i+1,q_star(i+1));
end
Path( (1:(T(n)-order)) + tacc ) = q_star;
tacc = tacc + T(n)-order;
end
end
end