spm_DFP.m 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727
  1. function [DEM] = spm_DFP(DEM)
  2. % Dynamic free-energy Fokker-Planck free-form scheme
  3. % FORMAT [DEM] = spm_DFP(DEM)
  4. %
  5. % DEM.M - hierarchical model
  6. % DEM.Y - output or data
  7. % DEM.U - inputs or prior expectation of causes
  8. % DEM.X - confounds
  9. %__________________________________________________________________________
  10. %
  11. % generative model
  12. %--------------------------------------------------------------------------
  13. % M(i).g = y(t) = g(x,v,P) {inline function, string or m-file}
  14. % M(i).f = dx/dt = f(x,v,P) {inline function, string or m-file}
  15. %
  16. % M(i).pE = prior expectation of p model-parameters
  17. % M(i).pC = prior covariances of p model-parameters
  18. % M(i).hE = prior expectation of h hyper-parameters (input noise)
  19. % M(i).hC = prior covariances of h hyper-parameters (input noise)
  20. % M(i).gE = prior expectation of g hyper-parameters (state noise)
  21. % M(i).gC = prior covariances of g hyper-parameters (state noise)
  22. % M(i).Q = precision components (input noise)
  23. % M(i).R = precision components (state noise)
  24. % M(i).V = fixed precision (input noise)
  25. % M(i).W = fixed precision (state noise)
  26. %
  27. % M(i).m = number of inputs v(i + 1);
  28. % M(i).n = number of states x(i);
  29. % M(i).l = number of output v(i);
  30. %
  31. % conditional moments of model-states - q(u)
  32. %--------------------------------------------------------------------------
  33. % qU.x = Conditional expectation of hidden states
  34. % qU.v = Conditional expectation of causal states
  35. % qU.e = Conditional residuals
  36. % qU.C = Conditional covariance: cov(v)
  37. % qU.S = Conditional covariance: cov(x)
  38. %
  39. % conditional moments of model-parameters - q(p)
  40. %--------------------------------------------------------------------------
  41. % qP.P = Conditional expectation
  42. % qP.Pi = Conditional expectation for each level
  43. % qP.C = Conditional covariance
  44. %
  45. % conditional moments of hyper-parameters (log-transformed) - q(h)
  46. %--------------------------------------------------------------------------
  47. % qH.h = Conditional expectation
  48. % qH.hi = Conditional expectation for each level
  49. % qH.C = Conditional covariance
  50. % qH.iC = Component precision: cov(vec(e[:})) = inv(kron(iC,iV))
  51. % qH.iV = Sequential precision
  52. %
  53. % F = log evidence = marginal likelihood = negative free energy
  54. %__________________________________________________________________________
  55. %
  56. % spm_DFP implements a variational Bayes (VB) scheme under the Laplace
  57. % approximation to the conditional densities of the model's, parameters (p)
  58. % and hyperparameters (h) of any analytic nonlinear hierarchical dynamic
  59. % model, with additive Gaussian innovations. It comprises three
  60. % variational steps (D,E and M) that update the conditional moments of u, p
  61. % and h respectively
  62. %
  63. % D: qu.u = max <L>q(p,h)
  64. % E: qp.p = max <L>q(u,h)
  65. % M: qh.h = max <L>q(u,p)
  66. %
  67. % where qu.u corresponds to the conditional expectation of hidden states x
  68. % and causal states v and so on. L is the ln p(y,u,p,h|M) under the model
  69. % M. The conditional covariances obtain analytically from the curvature of
  70. % L with respect to u, p and h.
  71. %
  72. % The D-step is implemented with variational filtering, which does not
  73. % assume a fixed form for the conditional density; it uses the sample
  74. % density of an ensemble of particles that drift up free-energy gradients
  75. % and 'explore' the local curvature though (Wiener) perturbations.
  76. %__________________________________________________________________________
  77. % Copyright (C) 2008 Wellcome Trust Centre for Neuroimaging
  78. % Karl Friston
  79. % $Id: spm_DFP.m 6540 2015-09-05 10:06:42Z karl $
  80. % Check model, data, priros and confounds and unpack
  81. %--------------------------------------------------------------------------
  82. [M,Y,U,X] = spm_DEM_set(DEM);
  83. MOVIE = 0;
  84. % find or create a DEM figure
  85. %--------------------------------------------------------------------------
  86. Fdem = spm_figure('GetWin','DEM');
  87. Fdfp = spm_figure('GetWin','DFP');
  88. % tolerance for changes in norm
  89. %--------------------------------------------------------------------------
  90. TOL = 1e-2;
  91. % order parameters (d = n = 1 for static models) and checks
  92. %==========================================================================
  93. d = M(1).E.d + 1; % embedding order of q(v)
  94. n = M(1).E.n + 1; % embedding order of q(x)
  95. s = M(1).E.s; % smoothness - s.d. of kernel (bins)
  96. try
  97. N = M(1).E.N; % number of particles
  98. catch
  99. N = 16; % number of particles
  100. end
  101. % number of states and parameters
  102. %--------------------------------------------------------------------------
  103. nY = size(Y,2); % number of samples
  104. nl = size(M,2); % number of levels
  105. ne = sum(cat(1,M.l)); % number of e (errors)
  106. nv = sum(cat(1,M.m)); % number of v (casual states)
  107. nx = sum(cat(1,M.n)); % number of x (hidden states)
  108. ny = M(1).l; % number of y (inputs)
  109. nc = M(end).l; % number of c (prior causes)
  110. nu = nv*d + nx*n; % number of generalised states
  111. kt = 1; % rate constant for D-Step
  112. % number of iterations
  113. %--------------------------------------------------------------------------
  114. if nx, nD = 1; else nD = 8; end
  115. try, nE = M(1).E.nE; catch, nE = 1; end
  116. try, nM = M(1).E.nM; catch, nM = 8; end
  117. try, nN = M(1).E.nN; catch, nN = 8; end
  118. % initialise regularisation parameters
  119. %--------------------------------------------------------------------------
  120. td = 1/nD; % integration time for D-Step
  121. te = 2; % integration time for E-Step
  122. % Precision (R) and covariance of generalised errors
  123. %--------------------------------------------------------------------------
  124. [iV,V] = spm_DEM_R(n,s);
  125. % precision components Q{} requiring [Re]ML estimators (M-Step)
  126. %==========================================================================
  127. Q = {};
  128. for i = 1:nl
  129. q0{i,i} = sparse(M(i).l,M(i).l);
  130. end
  131. for i = 1:nl - 1
  132. r0{i,i} = sparse(M(i).n,M(i).n);
  133. end
  134. Q0 = kron(iV,spm_cat(q0));
  135. R0 = kron(iV,spm_cat(r0));
  136. for i = 1:nl
  137. for j = 1:length(M(i).Q)
  138. q = q0;
  139. q{i,i} = M(i).Q{j};
  140. Q{end + 1} = blkdiag(kron(iV,spm_cat(q)),R0);
  141. end
  142. end
  143. for i = 1:nl - 1
  144. for j = 1:length(M(i).R)
  145. q = r0;
  146. q{i,i} = M(i).R{j};
  147. Q{end + 1} = blkdiag(Q0,kron(iV,spm_cat(q)));
  148. end
  149. end
  150. % and fixed components P
  151. %--------------------------------------------------------------------------
  152. Q0 = kron(iV,spm_cat(spm_diag({M.V})));
  153. R0 = kron(iV,spm_cat(spm_diag({M.W})));
  154. Qp = blkdiag(Q0,R0);
  155. nh = length(Q); % number of hyperparameters
  156. % hyperpriors
  157. %--------------------------------------------------------------------------
  158. ph.h = spm_vec({M.hE; M.gE}); % prior expectation of h
  159. ph.c = spm_cat(spm_diag({M.hC M.gC})); % prior covariances of h
  160. ph.ic = spm_pinv(ph.c); % prior precision
  161. qh.h = ph.h; % conditional expectation
  162. qh.c = ph.c; % conditional covariance
  163. % priors on parameters (in reduced parameter space)
  164. %==========================================================================
  165. pp.c = cell(nl,nl);
  166. qp.p = cell(nl,1);
  167. for i = 1:(nl - 1)
  168. % eigenvector reduction: p <- pE + qp.u*qp.p
  169. %----------------------------------------------------------------------
  170. qp.u{i} = spm_svd(M(i).pC); % basis for parameters
  171. M(i).p = size(qp.u{i},2); % number of qp.p
  172. qp.p{i} = sparse(M(i).p,1); % initial qp.p
  173. pp.c{i,i} = qp.u{i}'*M(i).pC*qp.u{i}; % prior covariance
  174. end
  175. Up = spm_cat(spm_diag(qp.u));
  176. % initialise and augment with confound parameters B; with flat priors
  177. %--------------------------------------------------------------------------
  178. np = sum(cat(1,M.p)); % number of model parameters
  179. nb = size(X,1); % number of confounds
  180. nn = nb*ny; % number of nuisance parameters
  181. nf = np + nn; % numer of free parameters
  182. ip = [1:np];
  183. ib = [1:nn] + np;
  184. pp.c = spm_cat(pp.c);
  185. pp.ic = spm_pinv(pp.c);
  186. % initialise conditional density q(p) (for D-Step)
  187. %--------------------------------------------------------------------------
  188. qp.e = spm_vec(qp.p);
  189. qp.c = sparse(nf,nf);
  190. qp.b = sparse(ny,nb);
  191. % initialise dedb
  192. %--------------------------------------------------------------------------
  193. for i = 1:nl
  194. dedbi{i,1} = sparse(M(i).l,nn);
  195. end
  196. for i = 1:nl - 1
  197. dndbi{i,1} = sparse(M(i).n,nn);
  198. end
  199. for i = 1:n
  200. dEdb{i,1} = spm_cat(dedbi);
  201. end
  202. for i = 1:n
  203. dNdb{i,1} = spm_cat(dndbi);
  204. end
  205. dEdb = [dEdb; dNdb];
  206. % initialise cell arrays for D-Step; e{i + 1} = (d/dt)^i[e] = e[i]
  207. %==========================================================================
  208. qu.x = cell(n + 1,1);
  209. qu.v = cell(n + 1,1);
  210. qy = cell(n + 1,1);
  211. qc = cell(n + 1,1);
  212. [qu.x{:}] = deal(sparse(nx,1));
  213. [qu.v{:}] = deal(sparse(nv,1));
  214. [qy{:} ] = deal(sparse(ny,1));
  215. [qc{:} ] = deal(sparse(nc,1));
  216. % initialise cell arrays for hierarchical structure of x[0] and v[0]
  217. %--------------------------------------------------------------------------
  218. x = {M(1:end - 1).x};
  219. v = {M(1 + 1:end).v};
  220. qu.x{1} = spm_vec(x);
  221. qu.v{1} = spm_vec(v);
  222. qu(1:N) = deal(qu);
  223. try xp = DEM.M(1).E.xp; catch, xp = 1; end
  224. try vp = DEM.M(1).E.vp; catch, vp = 1; end
  225. for i = 1:N
  226. qu(i).x{1} = qu(i).x{1} + randn(nx,1)/xp;
  227. qu(i).v{1} = qu(i).v{1} + randn(nv,1)/vp;
  228. end
  229. dq = {qu(1).x{1:n} qu(1).v{1:d} qy{1:n} qc{1:d}};
  230. % derivatives for Jacobian of D-step
  231. %--------------------------------------------------------------------------
  232. Dx = cell(n,n);
  233. Dv = cell(d,d);
  234. Dy = cell(n,n);
  235. Dc = cell(d,d);
  236. [Dx{:}] = deal(sparse(nx,nx));
  237. [Dv{:}] = deal(sparse(nv,nv));
  238. [Dy{:}] = deal(sparse(ny,ny));
  239. [Dc{:}] = deal(sparse(nc,nc));
  240. % Wiener process
  241. %--------------------------------------------------------------------------
  242. Ix = Dx;
  243. Iv = Dv;
  244. for i = 1:d
  245. Iv{i,i} = speye(nv,nv);
  246. end
  247. for i = 1:n
  248. Ix{i,i} = speye(nx,nx);
  249. end
  250. dfdw = spm_cat(spm_diag({Ix,Iv,Dy,Dc}));
  251. % add constant terms
  252. %--------------------------------------------------------------------------
  253. for i = 2:d
  254. Dv{i - 1,i} = speye(nv,nv);
  255. Dc{i - 1,i} = speye(nc,nc);
  256. end
  257. for i = 2:n
  258. Dx{i - 1,i} = speye(nx,nx);
  259. Dy{i - 1,i} = speye(ny,ny);
  260. end
  261. Du = spm_cat(spm_diag({Dx,Dv}));
  262. Dc = spm_cat(Dc);
  263. Dy = spm_cat(Dy);
  264. % gradients and curvatures for conditional uncertainty
  265. %--------------------------------------------------------------------------
  266. dUdu = sparse(nu,1);
  267. dUdp = sparse(nf,1);
  268. dUduu = sparse(nu,nu);
  269. dUdpp = sparse(nf,nf);
  270. % preclude unneceassry iterations
  271. %--------------------------------------------------------------------------
  272. if ~nh, nM = 1; end
  273. if ~nf, nE = 1; end
  274. if ~nf && ~nh, nN = 1; end
  275. % Iterate DEM
  276. %==========================================================================
  277. for iN = 1:nN
  278. % get time and celar persistent variables in evaluation routines
  279. %----------------------------------------------------------------------
  280. tic; clear spm_DEM_eval
  281. % E-Step: (with embedded D-Step)
  282. %======================================================================
  283. mp = zeros(nf,1);
  284. for iE = 1:nE
  285. % [re-]set accumulators for E-Step
  286. %------------------------------------------------------------------
  287. dFdp = zeros(nf,1);
  288. dFdpp = zeros(nf,nf);
  289. EE = sparse(0);
  290. ECE = sparse(0);
  291. qp.ic = sparse(0);
  292. qu_c = speye(1);
  293. % [re-]set precisions using ReML hyperparameter estimates
  294. %------------------------------------------------------------------
  295. iS = Qp;
  296. for i = 1:nh
  297. iS = iS + Q{i}*exp(qh.h(i));
  298. end
  299. % [re-]adjust for confounds
  300. %------------------------------------------------------------------
  301. Y = Y - qp.b*X;
  302. % [re-]set states & their derivatives
  303. %------------------------------------------------------------------
  304. try
  305. qu = QU{1};
  306. end
  307. % D-Step: (nD D-Steps for each sample)
  308. %==================================================================
  309. for iY = 1:nY
  310. % [re-]set states for static systems
  311. %--------------------------------------------------------------
  312. if ~nx
  313. try, qu = QU{iY}; end
  314. end
  315. % D-Step: until convergence for static systems
  316. %==============================================================
  317. for iD = 1:nD
  318. % sampling time
  319. %----------------------------------------------------------
  320. ts = iY + (iD - 1)/nD;
  321. % derivatives of responses and inputs
  322. %----------------------------------------------------------
  323. qy(1:n) = spm_DEM_embed(Y,n,ts);
  324. qc(1:d) = spm_DEM_embed(U,d,ts);
  325. % compute dEdb (derivatives of confounds)
  326. %----------------------------------------------------------
  327. b = spm_DEM_embed(X,n,ts);
  328. for i = 1:n
  329. dedbi{1} = -kron(b{i}',speye(ny,ny));
  330. dEdb{i,1} = spm_cat(dedbi);
  331. end
  332. % moments of ensemble density
  333. %==========================================================
  334. q = [qu.x];
  335. for i = 1:n + 1
  336. qx{i} = mean([q{i,:}],2);
  337. end
  338. q = [qu.v];
  339. for i = 1:n + 1
  340. qv{i} = mean([q{i,:}],2);
  341. end
  342. % mean field effects
  343. %----------------------------------------------------------
  344. dudt = spm_vec({qx(2:n + 1)
  345. qv(2:d + 1)
  346. qy(2:n + 1)
  347. qc(2:d + 1)});
  348. if iD == nD
  349. % ensemble covariance
  350. %------------------------------------------------------
  351. ux = [qu.x];
  352. uv = [qu.v];
  353. ux = ux(1:n,:);
  354. uv = uv(1:d,:);
  355. c = cov(spm_cat([ux; uv])');
  356. % quantities for M-Step
  357. %------------------------------------------------------
  358. quy.x = qx;
  359. quy.v = qv;
  360. quy.y = qy;
  361. quy.u = qc;
  362. [E,dE] = spm_DEM_eval(M,quy,qp);
  363. dE.dP = [dE.dp spm_cat(dEdb)];
  364. qu_c = qu_c*c;
  365. EE = E*E'+ EE;
  366. ECE = ECE + dE.du*c*dE.du'+ dE.dP*qp.c*dE.dP';
  367. % save states for qu(iY)
  368. %------------------------------------------------------
  369. qU(iY).x = qx;
  370. qU(iY).v = qv;
  371. qU(iY).y = qy;
  372. qU(iY).u = qc;
  373. qU(iY).e = E;
  374. qU(iY).c = c;
  375. end
  376. % evaluate functions:
  377. % e = v - g(x,v), dx/dt = f(x,v) and derivatives dE.dx, ...
  378. %==========================================================
  379. for iP = 1:N
  380. quy.x = qu(iP).x;
  381. quy.v = qu(iP).v;
  382. quy.y = qy;
  383. quy.u = qc;
  384. [e,de] = spm_DEM_eval(M,quy,qp);
  385. % conditional uncertainty about parameters
  386. %======================================================
  387. if np
  388. for i = 1:nu
  389. % 1st-order derivatives: dUdv, ... ;
  390. %----------------------------------------------
  391. CJ = qp.c(ip,ip)*de.dpu{i}'*iS;
  392. dUdu(i,1) = trace(CJ*de.dp);
  393. % 2nd-order derivatives
  394. %----------------------------------------------
  395. for j = 1:nu
  396. dUduu(i,j) = trace(CJ*de.dpu{j});
  397. end
  398. end
  399. end
  400. % D-step update: of causes v{i}, and other states u(i)
  401. %======================================================
  402. % compute dqdt: q = {u y c}; and dudt: u = {v{1:d} x}
  403. %------------------------------------------------------
  404. dIdu = -de.du'*iS*e - dUdu/2;
  405. % and second-order derivatives
  406. %------------------------------------------------------
  407. dIduu = -de.du'*iS*de.du - dUduu/2;
  408. dIduy = -de.du'*iS*de.dy;
  409. dIduc = -de.du'*iS*de.dc;
  410. % gradient
  411. %------------------------------------------------------
  412. dFdu = dudt;
  413. dFdu(1:nu) = dIdu + dFdu(1:nu);
  414. % Jacobian
  415. %------------------------------------------------------
  416. dFduu = spm_cat({dIduu dIduy dIduc;
  417. [] Dy [] ;
  418. [] [] Dc}) ;
  419. % update conditional modes of states
  420. %------------------------------------------------------
  421. du = spm_sde_dx(dFduu,dfdw,dFdu,td);
  422. dq = spm_unvec(du,dq);
  423. for i = 1:n
  424. qu(iP).x{i} = qu(iP).x{i} + dq{i};
  425. end
  426. for i = 1:d
  427. qu(iP).v{i} = qu(iP).v{i} + dq{i + n};
  428. end
  429. end
  430. end % D-Step
  431. % D-Step: save ensemble density and plot (over samples)
  432. %--------------------------------------------------------------
  433. QU{iY} = qu;
  434. figure(Fdfp)
  435. spm_DFP_plot(QU,nY)
  436. if MOVIE
  437. subplot(2,1,1)
  438. set(gca,'YLim',[-0.4 1.2])
  439. drawnow
  440. MOV(iY) = getframe(gca);
  441. end
  442. % Gradients and curvatures for E-Step:
  443. %==============================================================
  444. for i = ip
  445. % 1st-order derivatives: U = tr(C*J'*iS*J)
  446. %----------------------------------------------------------
  447. CJ = c*dE.dup{i}'*iS;
  448. dUdp(i,1) = trace(CJ*dE.du);
  449. % 2nd-order derivatives
  450. %----------------------------------------------------------
  451. for j = ip
  452. dUdpp(i,j) = trace(CJ*dE.dup{j});
  453. end
  454. end
  455. % Accumulate; dF/dP = <dL/dp>, dF/dpp = ...
  456. %--------------------------------------------------------------
  457. dFdp = dFdp - dUdp/2 - dE.dP'*iS*E;
  458. dFdpp = dFdpp - dUdpp/2 - dE.dP'*iS*dE.dP;
  459. qp.ic = qp.ic + dE.dP'*iS*dE.dP;
  460. end % sequence (iY)
  461. % augment with priors
  462. %------------------------------------------------------------------
  463. dFdp(ip) = dFdp(ip) - pp.ic*qp.e;
  464. dFdpp(ip,ip) = dFdpp(ip,ip) - pp.ic;
  465. qp.ic(ip,ip) = qp.ic(ip,ip) + pp.ic;
  466. qp.c = spm_pinv(qp.ic);
  467. % E-step: update expectation (p)
  468. %==================================================================
  469. % update conditional expectation
  470. %------------------------------------------------------------------
  471. dp = spm_dx(dFdpp,dFdp,{te});
  472. qp.e = qp.e + dp(ip);
  473. qp.p = spm_unvec(qp.e,qp.p);
  474. qp.b = spm_unvec(dp(ib),qp.b);
  475. mp = mp + dp;
  476. % convergence (E-Step)
  477. %------------------------------------------------------------------
  478. if (dFdp'*dp < 1e-2) || (norm(dp,1) < TOL), break, end
  479. end % E-Step
  480. % M-step - hyperparameters (h = exp(l))
  481. %======================================================================
  482. mh = zeros(nh,1);
  483. dFdh = zeros(nh,1);
  484. dFdhh = zeros(nh,nh);
  485. for iM = 1:nM
  486. % [re-]set precisions using ReML hyperparameter estimates
  487. %------------------------------------------------------------------
  488. iS = Qp;
  489. for i = 1:nh
  490. iS = iS + Q{i}*exp(qh.h(i));
  491. end
  492. S = inv(iS);
  493. dS = ECE + EE - S*nY;
  494. % 1st-order derivatives: dFdh = dF/dh
  495. %------------------------------------------------------------------
  496. for i = 1:nh
  497. dPdh{i} = Q{i}*exp(qh.h(i));
  498. dFdh(i,1) = -trace(dPdh{i}*dS)/2;
  499. end
  500. % 2nd-order derivatives: dFdhh
  501. %------------------------------------------------------------------
  502. for i = 1:nh
  503. for j = 1:nh
  504. dFdhh(i,j) = -trace(dPdh{i}*S*dPdh{j}*S*nY)/2;
  505. end
  506. end
  507. % hyperpriors
  508. %------------------------------------------------------------------
  509. qh.e = qh.h - ph.h;
  510. dFdh = dFdh - ph.ic*qh.e;
  511. dFdhh = dFdhh - ph.ic;
  512. % update ReML estimate of parameters
  513. %------------------------------------------------------------------
  514. dh = spm_dx(dFdhh,dFdh);
  515. qh.h = qh.h + dh;
  516. mh = mh + dh;
  517. % conditional covariance of hyperparameters
  518. %------------------------------------------------------------------
  519. qh.c = -spm_pinv(dFdhh);
  520. % convergence (M-Step)
  521. %------------------------------------------------------------------
  522. if (dFdh'*dh < 1e-2) || (norm(dh,1) < TOL), break, end
  523. end % M-Step
  524. % evaluate objective function (F)
  525. %======================================================================
  526. F(iN) = - trace(iS*EE)/2 ... % states (u)
  527. - trace(qp.e'*pp.ic*qp.e)/2 ... % parameters (p)
  528. - trace(qh.e'*ph.ic*qh.e)/2 ... % hyperparameters (h)
  529. + spm_logdet(qu_c)/2 ... % entropy q(u)
  530. + spm_logdet(qp.c)/2 ... % entropy q(p)
  531. + spm_logdet(qh.c)/2 ... % entropy q(h)
  532. - spm_logdet(pp.c)/2 ... % entropy - prior p
  533. - spm_logdet(ph.c)/2 ... % entropy - prior h
  534. + spm_logdet(iS)*nY/2 ... % entropy - error
  535. - n*ny*nY*log(2*pi)/2;
  536. % save model-states (for each time point)
  537. %======================================================================
  538. for t = 1:length(qU)
  539. v = spm_unvec(qU(t).v{1},v);
  540. x = spm_unvec(qU(t).x{1},x);
  541. e = spm_unvec(qU(t).e,{M.v});
  542. for i = 1:(nl - 1)
  543. Qu.v{i + 1}(:,t) = spm_vec(v{i});
  544. try
  545. Qu.x{i}(:,t) = spm_vec(x{i});
  546. end
  547. Qu.z{i}(:,t) = spm_vec(e{i});
  548. end
  549. Qu.v{1}(:,t) = spm_vec(qU(t).y{1} - e{1});
  550. Qu.z{nl}(:,t) = spm_vec(e{nl});
  551. % and conditional covariances
  552. %--------------------------------------------------------------
  553. i = [1:nx];
  554. Qu.S{t} = qU(t).c(i,i);
  555. i = [1:nv] + nx*n;
  556. Qu.C{t} = qU(t).c(i,i);
  557. end
  558. % report and break if convergence
  559. %------------------------------------------------------------------
  560. figure(Fdem)
  561. spm_DEM_qU(Qu)
  562. if np
  563. subplot(nl,4,4*nl)
  564. bar(full(Up*qp.e))
  565. xlabel({'parameters';'{minus prior}'})
  566. axis square, grid on
  567. end
  568. if length(F) > 2
  569. subplot(nl,4,4*nl - 1)
  570. plot(F(2:end))
  571. xlabel('iteractions')
  572. title('Log-evidence')
  573. axis square, grid on
  574. end
  575. drawnow
  576. % report (EM-Steps)
  577. %------------------------------------------------------------------
  578. str{1} = sprintf('DEM: %i (%i:%i:%i)',iN,iD,iE,iM);
  579. str{2} = sprintf('F:%.6e',full(F(iN)));
  580. str{3} = sprintf('p:%.2e',full(mp'*mp));
  581. str{4} = sprintf('h:%.2e',full(mh'*mh));
  582. str{5} = sprintf('(%.2e sec)',full(toc));
  583. fprintf('%-16s%-16s%-14s%-14s%-16s\n',str{:})
  584. if norm(mp) < TOL && norm(mh) < TOL, break, end
  585. end
  586. % Assemble output arguments
  587. %==========================================================================
  588. % conditional moments of model-parameters (rotated into original space)
  589. %--------------------------------------------------------------------------
  590. qP.P = spm_unvec(Up*qp.e + spm_vec(M.pE),M.pE);
  591. qP.C = Up*qp.c(ip,ip)*Up';
  592. qP.V = spm_unvec(diag(qP.C),M.pE);
  593. % conditional moments of hyper-parameters (log-transformed)
  594. %--------------------------------------------------------------------------
  595. qH.h = spm_unvec(qh.h,{{M.hE} {M.gE}});
  596. qH.g = qH.h{2};
  597. qH.h = qH.h{1};
  598. qH.C = qh.c;
  599. qH.V = spm_unvec(diag(qH.C),{{M.hE} {M.gE}});
  600. qH.W = qH.V{2};
  601. qH.V = qH.V{1};
  602. % assign output variables
  603. %--------------------------------------------------------------------------
  604. DEM.M = M;
  605. DEM.U = U; % causes
  606. DEM.X = X; % confounds
  607. DEM.QU = QU; % sample density of model-states
  608. DEM.qU = Qu; % conditional moments of model-states
  609. DEM.qP = qP; % conditional moments of model-parameters
  610. DEM.qH = qH; % conditional moments of hyper-parameters
  611. DEM.F = F; % [-ve] Free energy
  612. % set ButtonDownFcn
  613. %--------------------------------------------------------------------------
  614. if MOVIE
  615. figure(Fdfp), subplot(2,1,1)
  616. set(gca,'Userdata',{MOV,16})
  617. set(gca,'ButtonDownFcn','spm_DEM_ButtonDownFcn')
  618. end