spm_int_L.m 5.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169
  1. function [y] = spm_int_L(P,M,U,N)
  2. % Integrate a MIMO nonlinear system using a fixed Jacobian: J(x(0))
  3. % FORMAT [y] = spm_int_L(P,M,U,[N])
  4. % P - model parameters
  5. % M - model structure
  6. % U - input structure or matrix
  7. %
  8. % N - number of local linear iterations per time step [default: 1]
  9. %
  10. % y - (v x l) response y = g(x,u,P)
  11. %__________________________________________________________________________
  12. % Integrates the MIMO system described by
  13. %
  14. % dx/dt = f(x,u,P,M)
  15. % y = g(x,u,P,M)
  16. %
  17. % using the update scheme:
  18. %
  19. % x(t + dt) = x(t) + U*dx(t)/dt
  20. %
  21. % U = (expm(dt*J) - I)*inv(J)
  22. % J = df/dx
  23. %
  24. % at input times. This integration scheme evaluates the update matrix (U)
  25. % at the expansion point.
  26. %
  27. %--------------------------------------------------------------------------
  28. %
  29. % SPM solvers or integrators
  30. %
  31. % spm_int_ode: uses ode45 (or ode113) which are one and multi-step solvers
  32. % respectively. They can be used for any ODEs, where the Jacobian is
  33. % unknown or difficult to compute; however, they may be slow.
  34. %
  35. % spm_int_J: uses an explicit Jacobian-based update scheme that preserves
  36. % nonlinearities in the ODE: dx = (expm(dt*J) - I)*inv(J)*f. If the
  37. % equations of motion return J = df/dx, it will be used; otherwise it is
  38. % evaluated numerically, using spm_diff at each time point. This scheme is
  39. % infallible but potentially slow, if the Jacobian is not available (calls
  40. % spm_dx).
  41. %
  42. % spm_int_E: As for spm_int_J but uses the eigensystem of J(x(0)) to eschew
  43. % matrix exponentials and inversion during the integration. It is probably
  44. % the best compromise, if the Jacobian is not available explicitly.
  45. %
  46. % spm_int_B: As for spm_int_J but uses a first-order approximation to J
  47. % based on J(x(t)) = J(x(0)) + dJdx*x(t).
  48. %
  49. % spm_int_L: As for spm_int_B but uses J(x(0)).
  50. %
  51. % spm_int_U: like spm_int_J but only evaluates J when the input changes.
  52. % This can be useful if input changes are sparse (e.g., boxcar functions).
  53. % It is used primarily for integrating EEG models
  54. %
  55. % spm_int: Fast integrator that uses a bilinear approximation to the
  56. % Jacobian evaluated using spm_bireduce. This routine will also allow for
  57. % sparse sampling of the solution and delays in observing outputs. It is
  58. % used primarily for integrating fMRI models
  59. %__________________________________________________________________________
  60. % Copyright (C) 2008-2016 Wellcome Trust Centre for Neuroimaging
  61. % Karl Friston
  62. % $Id: spm_int_L.m 7143 2017-07-29 18:50:38Z karl $
  63. % convert U to U.u if necessary
  64. %--------------------------------------------------------------------------
  65. if ~isstruct(U), u.u = U; U = u; end
  66. try, dt = U.dt; catch, dt = 1; end
  67. if nargin < 4; N = 1; end
  68. % Initial states and inputs
  69. %--------------------------------------------------------------------------
  70. try
  71. x = M.x;
  72. catch
  73. x = sparse(0,1);
  74. M.x = x;
  75. end
  76. try
  77. u = U.u(1,:);
  78. catch
  79. u = sparse(1,M.m);
  80. end
  81. % add [0] states if not specified
  82. %--------------------------------------------------------------------------
  83. try
  84. f = spm_funcheck(M.f);
  85. catch
  86. f = @(x,u,P,M) sparse(0,1);
  87. M.n = 0;
  88. M.x = sparse(0,0);
  89. end
  90. M.f = f;
  91. % output nonlinearity, if specified
  92. %--------------------------------------------------------------------------
  93. try
  94. g = spm_funcheck(M.g);
  95. if isempty(g)
  96. g = @(x,u,P,M) x;
  97. end
  98. catch
  99. g = @(x,u,P,M) x;
  100. end
  101. M.g = g;
  102. % dx(t)/dt and Jacobian df/dx (and check for delay operator)
  103. %--------------------------------------------------------------------------
  104. D = 1;
  105. n = spm_length(x);
  106. if nargout(f) >= 3
  107. [fx,dfdx,D] = f(x,u,P,M);
  108. elseif nargout(f) == 2
  109. [fx,dfdx] = f(x,u,P,M);
  110. else
  111. dfdx = spm_cat(spm_diff(f,x,u,P,M,1));
  112. end
  113. % local linear update operator Q = (expm(dt*J) - I)*inv(J)
  114. %--------------------------------------------------------------------------
  115. dfdx = dfdx - speye(n,n)*exp(-16);
  116. Q = (spm_expm(dt*D*dfdx/N) - speye(n,n))/dfdx;
  117. % integrate
  118. %==========================================================================
  119. v = spm_vec(x);
  120. for i = 1:size(U.u,1)
  121. % input
  122. %----------------------------------------------------------------------
  123. u = U.u(i,:);
  124. try
  125. % update dx = (expm(dt*J) - I)*inv(J)*f(x,u)
  126. %------------------------------------------------------------------
  127. for j = 1:N
  128. v = v + Q*f(v,u,P,M);
  129. end
  130. % output - implement g(x)
  131. %------------------------------------------------------------------
  132. y(:,i) = g(v,u,P,M);
  133. catch
  134. % update dx = (expm(dt*J) - I)*inv(J)*f(x,u)
  135. %------------------------------------------------------------------
  136. for j = 1:N
  137. x = spm_vec(x) + Q*spm_vec(f(x,u,P,M));
  138. x = spm_unvec(x,M.x);
  139. end
  140. % output - implement g(x)
  141. %------------------------------------------------------------------
  142. y(:,i) = spm_vec(g(x,u,P,M));
  143. end
  144. end
  145. % transpose
  146. %--------------------------------------------------------------------------
  147. y = real(y');