spm_dcm_peb.m 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636
  1. function [PEB,P] = spm_dcm_peb(P,M,field)
  2. % Hierarchical (PEB) inversion of DCMs using BMR and VL
  3. % FORMAT [PEB,DCM] = spm_dcm_peb(DCM,M,field)
  4. % FORMAT [PEB,DCM] = spm_dcm_peb(DCM,X,field)
  5. %
  6. % DCM - {N [x M]} structure array of DCMs from N subjects
  7. % -------------------------------------------------------------------------
  8. % DCM{i}.M.pE - prior expectation of parameters
  9. % DCM{i}.M.pC - prior covariances of parameters
  10. % DCM{i}.Ep - posterior expectations
  11. % DCM{i}.Cp - posterior covariance
  12. % DCM{i}.F - free energy
  13. %
  14. % M.X - 2nd-level design matrix: X(:,1) = ones(N,1) [default]
  15. % M.bE - 3rd-level prior expectation [default: DCM{1}.M.pE]
  16. % M.bC - 3rd-level prior covariance [default: DCM{1}.M.pC/M.alpha]
  17. % M.pC - 2nd-level prior covariance [default: DCM{1}.M.pC/M.beta]
  18. % M.hE - 2nd-level prior expectation of log precisions [default: 0]
  19. % M.hC - 2nd-level prior covariances of log precisions [default: 1/16]
  20. % M.maxit - maximum iterations [default: 64]
  21. %
  22. % M.Q - covariance components: {'single','fields','all','none'}
  23. % M.alpha - optional scaling to specify M.bC [default = 1]
  24. % M.beta - optional scaling to specify M.pC [default = 16]
  25. % - if beta equals 0, sample variances will be used
  26. %
  27. % NB: the prior covariance of 2nd-level random effects is:
  28. % exp(M.hE)*DCM{1}.M.pC/M.beta [default DCM{1}.M.pC/16]
  29. %
  30. % NB2: to manually specify which parameters should be assigned to
  31. % which covariance components, M.Q can be set to a cell array of
  32. % [nxn] binary matrices, where n is the number of DCM
  33. % parameters. A value of M.Q{i}(n,n)==1 indicates that parameter
  34. % n should be modelled with component i.
  35. %
  36. % M.Xnames - cell array of names for second level parameters [default: {}]
  37. %
  38. % field - parameter fields in DCM{i}.Ep to optimise [default: {'A','B'}]
  39. % 'All' will invoke all fields. This argument effectively allows
  40. % one to specify the parameters that constitute random effects.
  41. %
  42. % PEB - hierarchical dynamic model
  43. % -------------------------------------------------------------------------
  44. % PEB.Snames - string array of first level model names
  45. % PEB.Pnames - string array of parameters of interest
  46. % PEB.Pind - indices of parameters in spm_vec(DCM{i}.Ep)
  47. % PEB.Xnames - names of second level parameters
  48. %
  49. % PEB.M.X - second level (between-subject) design matrix
  50. % PEB.M.W - second level (within-subject) design matrix
  51. % PEB.M.Q - precision [components] of second level random effects
  52. % PEB.M.pE - prior expectation of second level parameters
  53. % PEB.M.pC - prior covariance of second level parameters
  54. % PEB.M.hE - prior expectation of second level log-precisions
  55. % PEB.M.hC - prior covariance of second level log-precisions
  56. % PEB.Ep - posterior expectation of second level parameters
  57. % PEB.Eh - posterior expectation of second level log-precisions
  58. % PEB.Cp - posterior covariance of second level parameters
  59. % PEB.Ch - posterior covariance of second level log-precisions
  60. % PEB.Ce - expected covariance of second level random effects
  61. % PEB.F - free energy of second level model
  62. %
  63. % DCM - 1st level (reduced) DCM structures with empirical priors
  64. %
  65. % If DCM is an an (N x M} array, hierarchical inversion will be
  66. % applied to each model (i.e., each row) - and PEB will be a
  67. % {1 x M} cell array.
  68. %
  69. % This routine inverts a hierarchical DCM using variational Laplace and
  70. % Bayesian model reduction. In essence, it optimises the empirical priors
  71. % over the parameters of a set of first level DCMs, using second level or
  72. % between subject constraints specified in the design matrix X. This scheme
  73. % is efficient in the sense that it does not require inversion of the first
  74. % level DCMs - it just requires the prior and posterior densities from each
  75. % first level DCM to compute empirical priors under the implicit
  76. % hierarchical model. The output of this scheme (PEB) can be re-entered
  77. % recursively to invert deep hierarchical models. Furthermore, Bayesian
  78. % model comparison (BMC) can be specified in terms of the empirical priors
  79. % to perform BMC at the group level. Alternatively, subject-specific (first
  80. % level) posterior expectations can be used for classical inference in the
  81. % usual way. Note that these (summary statistics) are optimal in the sense
  82. % that they have been estimated under empirical (hierarchical) priors.
  83. %
  84. % If called with a single DCM, there are no between-subject effects and the
  85. % design matrix is assumed to model mixtures of parameters at the first
  86. % level. If called with a cell array, each column is assumed to contain 1st
  87. % level DCMs inverted under the same model.
  88. %
  89. %__________________________________________________________________________
  90. % Copyright (C) 2015-2016 Wellcome Trust Centre for Neuroimaging
  91. % Karl Friston
  92. % $Id: spm_dcm_peb.m 7476 2018-11-07 15:17:39Z peter $
  93. % get filenames and set up
  94. %==========================================================================
  95. if ~nargin
  96. [P, sts] = spm_select([2 Inf],'^DCM.*\.mat$','Select DCM*.mat files');
  97. if ~sts, return; end
  98. end
  99. if ischar(P), P = cellstr(P); end
  100. if isstruct(P), P = {P}; end
  101. % check for DEM structures
  102. %--------------------------------------------------------------------------
  103. try
  104. DEM = P;
  105. P = spm_dem2dcm(P);
  106. end
  107. % check parameter fields and design matrices
  108. %--------------------------------------------------------------------------
  109. try, load(P{1}); catch, DCM = P{1}; end
  110. if nargin < 2, M.X = ones(length(P),1); end
  111. if isnumeric(M), M = struct('X',M); end
  112. if ~isfield(M,'X'), M.X = ones(length(P),1); end
  113. if isempty(M.X), M.X = ones(length(P),1); end
  114. if nargin < 3; field = {'A','B'}; end
  115. if strcmpi(field,'all'); field = fieldnames(DCM.M.pE);end
  116. if ischar(field), field = {field}; end
  117. try, maxit = M.maxit; catch, maxit = 64; end
  118. % repeat for each model (column) if P is an array
  119. %==========================================================================
  120. if size(P,2) > 1
  121. for i = 1:size(P,2)
  122. [p,q] = spm_dcm_peb(P(:,i),M,field);
  123. PEB(i) = p;
  124. P(:,i) = q;
  125. end
  126. return
  127. end
  128. % get (first level) densities (summary statistics)
  129. %==========================================================================
  130. Ns = numel(P); % number of subjects
  131. if isfield(M,'bC') && Ns > 1
  132. q = spm_find_pC(M.bC,M.bE,field); % parameter indices
  133. elseif isnumeric(field)
  134. q = field; % parameter indices
  135. else
  136. q = spm_find_pC(DCM,field); % parameter indices
  137. end
  138. try
  139. Pstr = spm_fieldindices(DCM.M.pE,q); % field names
  140. catch
  141. if isfield(DCM,'Pnames')
  142. % PEB given as input. Field names have form covariate:fieldname
  143. %------------------------------------------------------------------
  144. Pstr = [];
  145. for i = 1:length(DCM.Xnames)
  146. str = strcat(DCM.Xnames{i}, ': ', DCM.Pnames);
  147. Pstr = [Pstr; str];
  148. end
  149. else
  150. % Generate field names
  151. %------------------------------------------------------------------
  152. Pstr = strcat('P', cellstr(num2str(q)));
  153. end
  154. end
  155. Np = numel(q); % number of parameters
  156. if Np == 1
  157. Pstr = {Pstr};
  158. end
  159. for i = 1:Ns
  160. % get first(within subject) level DCM
  161. %----------------------------------------------------------------------
  162. try, load(P{i}); catch, DCM = P{i}; end
  163. % posterior densities over all parameters
  164. %----------------------------------------------------------------------
  165. if isstruct(DCM.M.pC)
  166. pC{i} = diag(spm_vec(DCM.M.pC));
  167. else
  168. pC{i} = DCM.M.pC;
  169. end
  170. pC{i} = pC{i};
  171. pE{i} = spm_vec(DCM.M.pE);
  172. qE{i} = spm_vec(DCM.Ep);
  173. qC{i} = DCM.Cp;
  174. % deal with rank deficient priors
  175. %----------------------------------------------------------------------
  176. if i == 1
  177. PE = pE{i}(q);
  178. PC = pC{i}(q,q);
  179. U = spm_svd(PC);
  180. Ne = numel(pE{i});
  181. else
  182. if numel(pE{i}) ~= Ne
  183. error('Please ensure all DCMs have the same parameterisation');
  184. end
  185. end
  186. % select parameters in field
  187. %----------------------------------------------------------------------
  188. pE{i} = U'*pE{i}(q);
  189. pC{i} = U'*pC{i}(q,q)*U;
  190. qE{i} = U'*qE{i}(q);
  191. qC{i} = U'*qC{i}(q,q)*U;
  192. % shrink posterior to accommodate inefficient inversions
  193. %----------------------------------------------------------------------
  194. if Ns > 1
  195. qC{i} = spm_inv(spm_inv(qC{i}) + spm_inv(pC{i})/16);
  196. end
  197. % and free energy of model with full priors
  198. %----------------------------------------------------------------------
  199. iF(i) = DCM.F;
  200. end
  201. % hierarchical model design and defaults
  202. %==========================================================================
  203. % second level model
  204. %--------------------------------------------------------------------------
  205. if isfield(M,'alpha'), alpha = M.alpha; else, alpha = 1; end
  206. if isfield(M,'beta'), beta = M.beta; else, beta = 16; end
  207. if ~isfield(M,'W'), M.W = speye(Np,Np); end
  208. % covariance component specification
  209. %--------------------------------------------------------------------------
  210. Q = {};
  211. if ~isfield(M,'Q')
  212. OPTION = 'single';
  213. elseif iscell(M.Q) && isnumeric(M.Q{1})
  214. OPTION = 'manual';
  215. Q = M.Q;
  216. else
  217. OPTION = M.Q;
  218. end
  219. % design matrices
  220. %--------------------------------------------------------------------------
  221. if Ns > 1;
  222. % between-subject design matrices and prior expectations
  223. %======================================================================
  224. X = M.X;
  225. W = U'*M.W*U;
  226. else
  227. % within subject design
  228. %======================================================================
  229. OPTION = 'none';
  230. U = 1;
  231. X = 1;
  232. W = M.W;
  233. end
  234. % variable names
  235. %--------------------------------------------------------------------------
  236. if isfield(M,'Xnames')
  237. Xnames = M.Xnames;
  238. else
  239. Nx = size(X,2);
  240. Xnames = cell(1,Nx);
  241. for i = 1:Nx
  242. Xnames{i} = sprintf('Covariate %d',i);
  243. end
  244. end
  245. % get priors (from DCM if necessary) and ensure correct sizes
  246. %--------------------------------------------------------------------------
  247. if isfield(M,'bE')
  248. M.bE = spm_vec(M.bE);
  249. if size(M.bE,1) > Np && Ns > 1, M.bE = M.bE(q); end
  250. else
  251. M.bE = PE;
  252. end
  253. % prior covariance of (e.g.,) group effects
  254. %--------------------------------------------------------------------------
  255. if isfield(M,'bC')
  256. if isstruct(M.bC), M.bC = diag(spm_vec(M.bC)); end
  257. if size(M.bC,1) > Np && Ns > 1, M.bC = M.bC(q,q); end
  258. else
  259. M.bC = PC/alpha;
  260. end
  261. % between (e.g.,) subject covariances (for precision components Q)
  262. %--------------------------------------------------------------------------
  263. if isfield(M,'pC')
  264. if isstruct(M.pC), M.pC = diag(spm_vec(M.pC)); end
  265. if size(M.pC,1) > Np && Ns > 1, M.pC = M.pC(q,q); end
  266. elseif ~beta
  267. % If beta = 0, use variance of MAP estimators
  268. %----------------------------------------------------------------------
  269. M.pC = diag(var(spm_cat(qE),1,2));
  270. else
  271. % otherwise, use a scaled prior covariance
  272. %----------------------------------------------------------------------
  273. M.pC = PC/beta;
  274. end
  275. % prior precision (pP) and components (Q) for empirical covariance
  276. %--------------------------------------------------------------------------
  277. pQ = spm_inv(U'*M.pC*U);
  278. rP = pQ;
  279. switch OPTION
  280. case{'single'}
  281. % one between subject precision component
  282. %------------------------------------------------------------------
  283. Q = {pQ};
  284. case{'fields'}
  285. % between subject precision components (one for each field)
  286. %------------------------------------------------------------------
  287. pq = spm_inv(M.pC);
  288. for i = 1:length(field)
  289. j = spm_fieldindices(DCM.M.pE,field{i});
  290. j = find(ismember(q,j));
  291. Q{i} = sparse(Np,Np);
  292. Q{i}(j,j) = pq(j,j);
  293. Q{i} = U'*Q{i}*U;
  294. end
  295. case{'all'}
  296. % between subject precision components (one for each parameter)
  297. %------------------------------------------------------------------
  298. pq = spm_inv(M.pC);
  299. k = 1;
  300. for i = 1:Np
  301. qk = sparse(i,i,pq(i,i),Np,Np);
  302. qk = U'*qk*U;
  303. if any(qk(:))
  304. Q{k} = qk;
  305. k = k + 1;
  306. end
  307. end
  308. case {'manual'}
  309. % manually provided cell array of (binary) precision components
  310. %------------------------------------------------------------------
  311. pq = spm_inv(M.pC);
  312. k = 1;
  313. for i = 1:length(Q)
  314. j = find(diag(Q{i}));
  315. j = find(ismember(q,j));
  316. qk = sparse(Np,Np);
  317. qk(j,j) = pq(j,j);
  318. qk = U'*qk*U;
  319. if any(qk(:))
  320. Q{k} = qk;
  321. k = k + 1;
  322. end
  323. end
  324. case {'none'}
  325. % Do nothing
  326. otherwise
  327. warning('Unknown covariance component specification');
  328. end
  329. % number of parameters and effects
  330. %--------------------------------------------------------------------------
  331. Nx = size(X,2); % number of between subject effects
  332. Nw = size(W,2); % number of within subject effects
  333. Ng = numel(Q); % number of precision components
  334. Nb = Nw*Nx; % number of second level parameters
  335. % check for user-specified priors on log precision of second level effects
  336. %--------------------------------------------------------------------------
  337. gE = 0;
  338. gC = 1/16;
  339. try, gE = M.hE; end
  340. try, gC = M.hC; end
  341. try, Xc = M.bX; catch
  342. % adjust (second level) priors for the norm of explanatory variables
  343. %----------------------------------------------------------------------
  344. Xc = diag(size(X,1)./sum(X.^2));
  345. end
  346. % prior expectations and precisions of second level parameters
  347. %--------------------------------------------------------------------------
  348. Xe = spm_speye(Nx,1);
  349. bE = kron(Xe,U'*M.bE); % prior expectation of group effects
  350. gE = zeros(Ng,1) + gE; % prior expectation of log precision
  351. bC = kron(Xc,U'*M.bC*U); % prior covariance of group effects
  352. gC = eye(Ng,Ng)*gC; % prior covariance of log precision
  353. bP = spm_inv(bC);
  354. gP = spm_inv(gC);
  355. % initialise parameters
  356. %--------------------------------------------------------------------------
  357. b = bE;
  358. g = gE;
  359. ipC = spm_cat({bP [];
  360. [] gP});
  361. % variational Laplace
  362. %--------------------------------------------------------------------------
  363. t = -4; % Fisher scoring parameter
  364. for n = 1:maxit
  365. % compute prior precision (with a lower bound of pQ/exp(8))
  366. %----------------------------------------------------------------------
  367. if Ng > 0
  368. rP = pQ*exp(-8);
  369. for i = 1:Ng
  370. rP = rP + exp(g(i))*Q{i};
  371. end
  372. end
  373. rC = spm_inv(rP);
  374. % update model parameters
  375. %======================================================================
  376. % Gradient and curvature with respect to free energy
  377. %----------------------------------------------------------------------
  378. F = 0;
  379. dFdb = -bP*(b - bE);
  380. dFdbb = -bP;
  381. dFdg = -gP*(g - gE);
  382. dFdgg = -gP;
  383. dFdbg = zeros(Nb,Ng);
  384. for i = 1:Ns
  385. % get empirical prior expectations and reduced 1st level posterior
  386. %------------------------------------------------------------------
  387. Xi = kron(X(i,:),W);
  388. rE{i} = Xi*b;
  389. [Fi,sE,sC] = spm_log_evidence_reduce(qE{i},qC{i},pE{i},pC{i},rE{i},rC);
  390. % supplement gradients and curvatures
  391. %------------------------------------------------------------------
  392. F = F + Fi + iF(i);
  393. dE = sE - rE{i};
  394. dFdb = dFdb + Xi'*rP*dE;
  395. dFdbb = dFdbb + Xi'*(rP*sC*rP - rP)*Xi;
  396. for j = 1:Ng
  397. dFdgj = exp(g(j))*(trace((rC - sC)*Q{j}) - dE'*Q{j}*dE)/2;
  398. dFdg(j) = dFdg(j) + dFdgj;
  399. dFdgg(j,j) = dFdgg(j,j) + dFdgj;
  400. dFdbgj = exp(g(j))*(Xi - sC*rP*Xi)'*Q{j}*dE;
  401. dFdbg(:,j) = dFdbg(:,j) + dFdbgj;
  402. for k = 1:Ng
  403. dFdggj = exp(g(j) + g(k))*(trace((rC*Q{k}*rC - sC*Q{k}*sC)*Q{j})/2 - dE'*Q{k}*sC*Q{j}*dE);
  404. dFdgg(j,k) = dFdgg(j,k) - dFdggj;
  405. end
  406. end
  407. end
  408. % Free-energy and update parameters
  409. %======================================================================
  410. dFdp = [dFdb; dFdg];
  411. dFdpp = spm_cat({dFdbb dFdbg;
  412. dFdbg' dFdgg});
  413. Cp = spm_inv(-dFdpp);
  414. Cb = spm_inv(-dFdbb);
  415. Cg = spm_inv(-dFdgg);
  416. % second level complexity component of free energy
  417. %----------------------------------------------------------------------
  418. Fb = b'*bP*b;
  419. Fg = g'*gP*g;
  420. Fc = Fb/2 + Fg/2 - spm_logdet(ipC*Cp)/2;
  421. F = F - Fc;
  422. % best free energy so far
  423. %----------------------------------------------------------------------
  424. if n == 1, F0 = F; end
  425. % convergence and regularisation
  426. %======================================================================
  427. % if F is increasing save current expansion point
  428. %----------------------------------------------------------------------
  429. if F >= F0
  430. dF = F - F0;
  431. F0 = F;
  432. save('tmp.mat','b','g','F0','dFdb','dFdbb','dFdg','dFdgg');
  433. % decrease regularisation
  434. %------------------------------------------------------------------
  435. t = min(t + 1/4,2);
  436. else
  437. % otherwise, retrieve expansion point and increase regularisation
  438. %------------------------------------------------------------------
  439. t = max(t - 1,-4);
  440. load('tmp.mat');
  441. end
  442. % Fisher scoring
  443. %----------------------------------------------------------------------
  444. dp = spm_dx(dFdpp,dFdp,{t});
  445. % Suppress conditional dependencies for stable convergence
  446. %------------------------------------------------------------------
  447. if norm(dp,1) < 8, else
  448. dFdpp = spm_cat({dFdbb [];
  449. [] dFdgg});
  450. dp = spm_dx(dFdpp,dFdp,{t});
  451. end
  452. [db,dg] = spm_unvec(dp,b,g);
  453. b = b + db;
  454. g = g + tanh(dg);
  455. % Convergence
  456. %======================================================================
  457. if ~isfield(M,'noplot')
  458. fprintf('VL Iteration %-8d: F = %-3.2f dF: %2.4f [%+2.2f]\n',n,full(F),full(dF),t);
  459. end
  460. if (n > 4) && (t <= -4 || dF < 1e-4)
  461. fprintf('VL Iteration %-8d: F = %-3.2f dF: %2.4f [%+2.2f]\n',n,full(F),full(dF),t);
  462. break
  463. end
  464. end
  465. % assemble output structure
  466. %==========================================================================
  467. for i = 1:Ns
  468. % get first(within subject) level DCM
  469. %----------------------------------------------------------------------
  470. try
  471. load(P{i},'DCM');
  472. Sstr{i} = P{i};
  473. catch
  474. DCM = P{i};
  475. try
  476. Sstr{i} = DCM.name;
  477. catch
  478. Sstr{i} = sprintf('Subject %i',i);
  479. end
  480. end
  481. % evaluate reduced first level parameters if required
  482. %----------------------------------------------------------------------
  483. if nargout > 1
  484. % posterior densities over all parameters
  485. %------------------------------------------------------------------
  486. if isstruct(DCM.M.pC)
  487. pC{i} = diag(spm_vec(DCM.M.pC));
  488. else
  489. pC{i} = DCM.M.pC;
  490. end
  491. pE{i} = DCM.M.pE;
  492. qE{i} = DCM.Ep;
  493. qC{i} = DCM.Cp;
  494. % augment empirical priors
  495. %------------------------------------------------------------------
  496. RP = spm_inv(pC{i});
  497. RP(q,q) = U*rP*U';
  498. RC = spm_inv(RP);
  499. RE = spm_vec(pE{i});
  500. RE(q) = U*rE{i};
  501. RE = spm_unvec(RE,pE{i});
  502. % First level BMR (supplemented with second level complexity)
  503. %------------------------------------------------------------------
  504. [Fi,sE,sC] = spm_log_evidence_reduce(qE{i},qC{i},pE{i},pC{i},RE,RC);
  505. DCM.M.pE = RE;
  506. DCM.M.pC = RC;
  507. DCM.Ep = sE;
  508. DCM.Cp = sC;
  509. DCM.F = Fi + iF(i) - Fc/Ns;
  510. P{i} = DCM;
  511. end
  512. end
  513. % second level results
  514. %--------------------------------------------------------------------------
  515. PEB.Snames = Sstr';
  516. PEB.Pnames = Pstr';
  517. PEB.Pind = q;
  518. PEB.Xnames = Xnames;
  519. Ub = kron(eye(Nx,Nx),U);
  520. PEB.M.X = X;
  521. PEB.M.W = W;
  522. PEB.M.U = U;
  523. PEB.M.pE = kron(Xe,M.bE);
  524. PEB.M.pC = kron(Xc,M.bC);
  525. PEB.M.hE = gE;
  526. PEB.M.hC = gC;
  527. PEB.Ep = U*reshape(b,Nw,Nx);
  528. PEB.Eh = g;
  529. PEB.Ch = Cg;
  530. PEB.Cp = Ub*Cb*Ub';
  531. PEB.Ce = U*rC*U';
  532. PEB.F = F;
  533. for i = 1:length(Q)
  534. PEB.M.Q{i} = U*Q{i}*U';
  535. end
  536. spm_unlink('tmp.mat');
  537. % check for DEM structures
  538. %--------------------------------------------------------------------------
  539. try, P = spm_dem2dcm(DEM,P); end