current_source_density.py 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332
  1. # -*- coding: utf-8 -*-
  2. """'Current Source Density analysis (CSD) is a class of methods of analysis of
  3. extracellular electric potentials recorded at multiple sites leading to
  4. estimates of current sources generating the measured potentials. It is usually
  5. applied to low-frequency part of the potential (called the Local Field
  6. Potential, LFP) and to simultaneous recordings or to recordings taken with
  7. fixed time reference to the onset of specific stimulus (Evoked Potentials)'
  8. (Definition by Prof.Daniel K. Wójcik for Encyclopedia of Computational
  9. Neuroscience)
  10. CSD is also called as Source Localization or Source Imaging in the EEG circles.
  11. Here are CSD methods for different types of electrode configurations.
  12. 1D - laminar probe like electrodes.
  13. 2D - Microelectrode Array like
  14. 3D - UtahArray or multiple laminar probes.
  15. The following methods have been implemented so far
  16. 1D - StandardCSD, DeltaiCSD, SplineiCSD, StepiCSD, KCSD1D
  17. 2D - KCSD2D, MoIKCSD (Saline layer on top of slice)
  18. 3D - KCSD3D
  19. Each of these methods listed have some advantages. The KCSD methods for
  20. instance can handle broken or irregular electrode configurations electrode
  21. Keywords: LFP; CSD; Multielectrode; Laminar electrode; Barrel cortex
  22. Citation Policy: See ./current_source_density_src/README.md
  23. Contributors to this current source density estimation module are:
  24. Chaitanya Chintaluri(CC), Espen Hagen(EH) and Michał Czerwinski(MC).
  25. EH implemented the iCSD methods and StandardCSD
  26. CC implemented the kCSD methods, kCSD1D(MC and CC)
  27. CC and EH developed the interface to elephant.
  28. """
  29. from __future__ import division
  30. import neo
  31. import quantities as pq
  32. import numpy as np
  33. from scipy import io
  34. from scipy.integrate import simps
  35. from elephant.current_source_density_src import KCSD
  36. from elephant.current_source_density_src import icsd
  37. import elephant.current_source_density_src.utility_functions as utils
  38. utils.patch_quantities()
  39. available_1d = ['StandardCSD', 'DeltaiCSD', 'StepiCSD', 'SplineiCSD', 'KCSD1D']
  40. available_2d = ['KCSD2D', 'MoIKCSD']
  41. available_3d = ['KCSD3D']
  42. kernel_methods = ['KCSD1D', 'KCSD2D', 'KCSD3D', 'MoIKCSD']
  43. icsd_methods = ['DeltaiCSD', 'StepiCSD', 'SplineiCSD']
  44. py_iCSD_toolbox = ['StandardCSD'] + icsd_methods
  45. def estimate_csd(lfp, coords=None, method=None,
  46. process_estimate=True, **kwargs):
  47. """
  48. Fuction call to compute the current source density (CSD) from extracellular
  49. potential recordings(local-field potentials - LFP) using laminar electrodes
  50. or multi-contact electrodes with 2D or 3D geometries.
  51. Parameters
  52. ----------
  53. lfp : list(neo.AnalogSignal type objects)
  54. positions of electrodes can be added as neo.RecordingChannel
  55. coordinate or sent externally as a func argument (See coords)
  56. coords : [Optional] corresponding spatial coordinates of the electrodes
  57. Defaults to None
  58. Otherwise looks for RecordingChannels coordinate
  59. method : string
  60. Pick a method corresonding to the setup, in this implementation
  61. For Laminar probe style (1D), use 'KCSD1D' or 'StandardCSD',
  62. or 'DeltaiCSD' or 'StepiCSD' or 'SplineiCSD'
  63. For MEA probe style (2D), use 'KCSD2D', or 'MoIKCSD'
  64. For array of laminar probes (3D), use 'KCSD3D'
  65. Defaults to None
  66. process_estimate : bool
  67. In the py_iCSD_toolbox this corresponds to the filter_csd -
  68. the parameters are passed as kwargs here ie., f_type and f_order
  69. In the kcsd methods this corresponds to cross_validate -
  70. the parameters are passed as kwargs here ie., lambdas and Rs
  71. Defaults to True
  72. kwargs : parameters to each method
  73. The parameters corresponding to the method chosen
  74. See the documentation of the individual method
  75. Default is {} - picks the best parameters,
  76. Returns
  77. -------
  78. Estimated CSD
  79. neo.AnalogSignal Object
  80. annotated with the spatial coordinates
  81. Raises
  82. ------
  83. AttributeError
  84. No units specified for electrode spatial coordinates
  85. ValueError
  86. Invalid function arguments, wrong method name, or
  87. mismatching coordinates
  88. TypeError
  89. Invalid cv_param argument passed
  90. """
  91. if not isinstance(lfp, neo.AnalogSignal):
  92. raise TypeError('Parameter `lfp` must be a list(neo.AnalogSignal \
  93. type objects')
  94. if coords is None:
  95. coords = lfp.channel_index.coordinates
  96. # for ii in lfp:
  97. # coords.append(ii.channel_index.coordinate.rescale(pq.mm))
  98. else:
  99. scaled_coords = []
  100. for coord in coords:
  101. try:
  102. scaled_coords.append(coord.rescale(pq.mm))
  103. except AttributeError:
  104. raise AttributeError('No units given for electrode spatial \
  105. coordinates')
  106. coords = scaled_coords
  107. if method is None:
  108. raise ValueError('Must specify a method of CSD implementation')
  109. if len(coords) != len(lfp):
  110. raise ValueError('Number of signals and coords is not same')
  111. for ii in coords: # CHECK for Dimensionality of electrodes
  112. if len(ii) > 3:
  113. raise ValueError('Invalid number of coordinate positions')
  114. dim = len(coords[0]) # TODO : Generic co-ordinates!
  115. if dim == 1 and (method not in available_1d):
  116. raise ValueError('Invalid method, Available options are:',
  117. available_1d)
  118. if dim == 2 and (method not in available_2d):
  119. raise ValueError('Invalid method, Available options are:',
  120. available_2d)
  121. if dim == 3 and (method not in available_3d):
  122. raise ValueError('Invalid method, Available options are:',
  123. available_3d)
  124. if method in kernel_methods:
  125. input_array = np.zeros((len(lfp), lfp[0].magnitude.shape[0]))
  126. for ii, jj in enumerate(lfp):
  127. input_array[ii, :] = jj.rescale(pq.mV).magnitude
  128. kernel_method = getattr(KCSD, method) # fetch the class 'KCSD1D'
  129. lambdas = kwargs.pop('lambdas', None)
  130. Rs = kwargs.pop('Rs', None)
  131. k = kernel_method(np.array(coords), input_array, **kwargs)
  132. if process_estimate:
  133. k.cross_validate(lambdas, Rs)
  134. estm_csd = k.values()
  135. estm_csd = np.rollaxis(estm_csd, -1, 0)
  136. output = neo.AnalogSignal(estm_csd * pq.uA / pq.mm**3,
  137. t_start=lfp.t_start,
  138. sampling_rate=lfp.sampling_rate)
  139. if dim == 1:
  140. output.annotate(x_coords=k.estm_x)
  141. elif dim == 2:
  142. output.annotate(x_coords=k.estm_x, y_coords=k.estm_y)
  143. elif dim == 3:
  144. output.annotate(x_coords=k.estm_x, y_coords=k.estm_y,
  145. z_coords=k.estm_z)
  146. elif method in py_iCSD_toolbox:
  147. coords = np.array(coords) * coords[0].units
  148. if method in icsd_methods:
  149. try:
  150. coords = coords.rescale(kwargs['diam'].units)
  151. except KeyError: # Then why specify as a default in icsd?
  152. # All iCSD methods explicitly assume a source
  153. # diameter in contrast to the stdCSD that
  154. # implicitly assume infinite source radius
  155. raise ValueError("Parameter diam must be specified for iCSD \
  156. methods: {}".format(", ".join(icsd_methods)))
  157. if 'f_type' in kwargs:
  158. if (kwargs['f_type'] is not 'identity') and \
  159. (kwargs['f_order'] is None):
  160. raise ValueError("The order of {} filter must be \
  161. specified".format(kwargs['f_type']))
  162. lfp = neo.AnalogSignal(np.asarray(lfp).T, units=lfp.units,
  163. sampling_rate=lfp.sampling_rate)
  164. csd_method = getattr(icsd, method) # fetch class from icsd.py file
  165. csd_estimator = csd_method(lfp=lfp.magnitude.T * lfp.units,
  166. coord_electrode=coords.flatten(),
  167. **kwargs)
  168. csd_pqarr = csd_estimator.get_csd()
  169. if process_estimate:
  170. csd_pqarr_filtered = csd_estimator.filter_csd(csd_pqarr)
  171. output = neo.AnalogSignal(csd_pqarr_filtered.T,
  172. t_start=lfp.t_start,
  173. sampling_rate=lfp.sampling_rate)
  174. else:
  175. output = neo.AnalogSignal(csd_pqarr.T, t_start=lfp.t_start,
  176. sampling_rate=lfp.sampling_rate)
  177. output.annotate(x_coords=coords)
  178. return output
  179. def generate_lfp(csd_profile, ele_xx, ele_yy=None, ele_zz=None,
  180. xlims=[0., 1.], ylims=[0., 1.], zlims=[0., 1.], res=50):
  181. """Forward modelling for the getting the potentials for testing CSD
  182. Parameters
  183. ----------
  184. csd_profile : fuction that computes True CSD profile
  185. Available options are (see ./csd/utility_functions.py)
  186. 1D : gauss_1d_dipole
  187. 2D : large_source_2D and small_source_2D
  188. 3D : gauss_3d_dipole
  189. ele_xx : np.array
  190. Positions of the x coordinates of the electrodes
  191. ele_yy : np.array
  192. Positions of the y coordinates of the electrodes
  193. Defaults ot None, use in 2D or 3D cases only
  194. ele_zz : np.array
  195. Positions of the z coordinates of the electrodes
  196. Defaults ot None, use in 3D case only
  197. x_lims : [start, end]
  198. The starting spatial coordinate and the ending for integration
  199. Defaults to [0.,1.]
  200. y_lims : [start, end]
  201. The starting spatial coordinate and the ending for integration
  202. Defaults to [0.,1.], use only in 2D and 3D case
  203. z_lims : [start, end]
  204. The starting spatial coordinate and the ending for integration
  205. Defaults to [0.,1.], use only in 3D case
  206. res : int
  207. The resolution of the integration
  208. Defaults to 50
  209. Returns
  210. -------
  211. LFP : list(neo.AnalogSignal type objects)
  212. The potentials created by the csd profile at the electrode positions
  213. The electrode postions are attached as RecordingChannel's coordinate
  214. """
  215. def integrate_1D(x0, csd_x, csd, h):
  216. m = np.sqrt((csd_x - x0)**2 + h**2) - abs(csd_x - x0)
  217. y = csd * m
  218. I = simps(y, csd_x)
  219. return I
  220. def integrate_2D(x, y, xlin, ylin, csd, h, X, Y):
  221. Ny = ylin.shape[0]
  222. m = np.sqrt((x - X)**2 + (y - Y)**2)
  223. m[m < 0.0000001] = 0.0000001
  224. y = np.arcsinh(2 * h / m) * csd
  225. I = np.zeros(Ny)
  226. for i in range(Ny):
  227. I[i] = simps(y[:, i], ylin)
  228. F = simps(I, xlin)
  229. return F
  230. def integrate_3D(x, y, z, xlim, ylim, zlim, csd, xlin, ylin, zlin,
  231. X, Y, Z):
  232. Nz = zlin.shape[0]
  233. Ny = ylin.shape[0]
  234. m = np.sqrt((x - X)**2 + (y - Y)**2 + (z - Z)**2)
  235. m[m < 0.0000001] = 0.0000001
  236. z = csd / m
  237. Iy = np.zeros(Ny)
  238. for j in range(Ny):
  239. Iz = np.zeros(Nz)
  240. for i in range(Nz):
  241. Iz[i] = simps(z[:, j, i], zlin)
  242. Iy[j] = simps(Iz, ylin)
  243. F = simps(Iy, xlin)
  244. return F
  245. dim = 1
  246. if ele_zz is not None:
  247. dim = 3
  248. elif ele_yy is not None:
  249. dim = 2
  250. x = np.linspace(xlims[0], xlims[1], res)
  251. if dim >= 2:
  252. y = np.linspace(ylims[0], ylims[1], res)
  253. if dim == 3:
  254. z = np.linspace(zlims[0], zlims[1], res)
  255. sigma = 1.0
  256. h = 50.
  257. pots = np.zeros(len(ele_xx))
  258. if dim == 1:
  259. chrg_x = np.linspace(xlims[0], xlims[1], res)
  260. csd = csd_profile(chrg_x)
  261. for ii in range(len(ele_xx)):
  262. pots[ii] = integrate_1D(ele_xx[ii], chrg_x, csd, h)
  263. pots /= 2. * sigma # eq.: 26 from Potworowski et al
  264. ele_pos = ele_xx
  265. elif dim == 2:
  266. chrg_x, chrg_y = np.mgrid[xlims[0]:xlims[1]:np.complex(0, res),
  267. ylims[0]:ylims[1]:np.complex(0, res)]
  268. csd = csd_profile(chrg_x, chrg_y)
  269. for ii in range(len(ele_xx)):
  270. pots[ii] = integrate_2D(ele_xx[ii], ele_yy[ii],
  271. x, y, csd, h, chrg_x, chrg_y)
  272. pots /= 2 * np.pi * sigma
  273. ele_pos = np.vstack((ele_xx, ele_yy)).T
  274. elif dim == 3:
  275. chrg_x, chrg_y, chrg_z = np.mgrid[xlims[0]:xlims[1]:np.complex(0, res),
  276. ylims[0]:ylims[1]:np.complex(0, res),
  277. zlims[0]:zlims[1]:np.complex(0, res)]
  278. csd = csd_profile(chrg_x, chrg_y, chrg_z)
  279. xlin = chrg_x[:, 0, 0]
  280. ylin = chrg_y[0, :, 0]
  281. zlin = chrg_z[0, 0, :]
  282. for ii in range(len(ele_xx)):
  283. pots[ii] = integrate_3D(ele_xx[ii], ele_yy[ii], ele_zz[ii],
  284. xlims, ylims, zlims, csd,
  285. xlin, ylin, zlin,
  286. chrg_x, chrg_y, chrg_z)
  287. pots /= 4 * np.pi * sigma
  288. ele_pos = np.vstack((ele_xx, ele_yy, ele_zz)).T
  289. pots = np.reshape(pots, (-1, 1)) * pq.mV
  290. ele_pos = ele_pos * pq.mm
  291. lfp = []
  292. ch = neo.ChannelIndex(index=range(len(pots)))
  293. for ii in range(len(pots)):
  294. lfp.append(pots[ii])
  295. # lfp = neo.AnalogSignal(lfp, sampling_rate=1000*pq.Hz, units='mV')
  296. asig = neo.AnalogSignal(lfp, sampling_rate=pq.kHz, units='mV')
  297. ch.coordinates = ele_pos
  298. ch.analogsignals.append(asig)
  299. ch.create_relationship()
  300. return asig