current_source_density.py 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340
  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, print_function, unicode_literals
  30. import neo
  31. import numpy as np
  32. import quantities as pq
  33. from scipy.integrate import simps
  34. import elephant.current_source_density_src.utility_functions as utils
  35. from elephant.current_source_density_src import KCSD, icsd
  36. from elephant.utils import deprecated_alias
  37. __all__ = [
  38. "estimate_csd",
  39. "generate_lfp"
  40. ]
  41. utils.patch_quantities()
  42. available_1d = ['StandardCSD', 'DeltaiCSD', 'StepiCSD', 'SplineiCSD', 'KCSD1D']
  43. available_2d = ['KCSD2D', 'MoIKCSD']
  44. available_3d = ['KCSD3D']
  45. kernel_methods = ['KCSD1D', 'KCSD2D', 'KCSD3D', 'MoIKCSD']
  46. icsd_methods = ['DeltaiCSD', 'StepiCSD', 'SplineiCSD']
  47. py_iCSD_toolbox = ['StandardCSD'] + icsd_methods
  48. @deprecated_alias(coords='coordinates')
  49. def estimate_csd(lfp, coordinates=None, method=None,
  50. process_estimate=True, **kwargs):
  51. """
  52. Function call to compute the current source density (CSD) from
  53. extracellular potential recordings(local-field potentials - LFP) using
  54. laminar electrodes or multi-contact electrodes with 2D or 3D geometries.
  55. Parameters
  56. ----------
  57. lfp : neo.AnalogSignal
  58. positions of electrodes can be added as neo.RecordingChannel
  59. coordinate or sent externally as a func argument (See coords)
  60. coordinates : [Optional] corresponding spatial coordinates of the
  61. electrodes.
  62. Defaults to None
  63. Otherwise looks for ChannelIndex coordinate
  64. method : string
  65. Pick a method corresponding to the setup, in this implementation
  66. For Laminar probe style (1D), use 'KCSD1D' or 'StandardCSD',
  67. or 'DeltaiCSD' or 'StepiCSD' or 'SplineiCSD'
  68. For MEA probe style (2D), use 'KCSD2D', or 'MoIKCSD'
  69. For array of laminar probes (3D), use 'KCSD3D'
  70. Defaults to None
  71. process_estimate : bool
  72. In the py_iCSD_toolbox this corresponds to the filter_csd -
  73. the parameters are passed as kwargs here ie., f_type and f_order
  74. In the kcsd methods this corresponds to cross_validate -
  75. the parameters are passed as kwargs here ie., lambdas and Rs
  76. Defaults to True
  77. kwargs : parameters to each method
  78. The parameters corresponding to the method chosen
  79. See the documentation of the individual method
  80. Default is {} - picks the best parameters,
  81. Returns
  82. -------
  83. Estimated CSD
  84. neo.AnalogSignal object
  85. annotated with the spatial coordinates
  86. Raises
  87. ------
  88. AttributeError
  89. No units specified for electrode spatial coordinates
  90. ValueError
  91. Invalid function arguments, wrong method name, or
  92. mismatching coordinates
  93. TypeError
  94. Invalid cv_param argument passed
  95. """
  96. if not isinstance(lfp, neo.AnalogSignal):
  97. raise TypeError('Parameter `lfp` must be a neo.AnalogSignal object')
  98. if coordinates is None:
  99. coordinates = lfp.channel_index.coordinates
  100. else:
  101. scaled_coords = []
  102. for coord in coordinates:
  103. try:
  104. scaled_coords.append(coord.rescale(pq.mm))
  105. except AttributeError:
  106. raise AttributeError('No units given for electrode spatial \
  107. coordinates')
  108. coordinates = scaled_coords
  109. if method is None:
  110. raise ValueError('Must specify a method of CSD implementation')
  111. if len(coordinates) != lfp.shape[1]:
  112. raise ValueError('Number of signals and coords is not same')
  113. for ii in coordinates: # CHECK for Dimensionality of electrodes
  114. if len(ii) > 3:
  115. raise ValueError('Invalid number of coordinate positions')
  116. dim = len(coordinates[0]) # TODO : Generic co-ordinates!
  117. if dim == 1 and (method not in available_1d):
  118. raise ValueError('Invalid method, Available options are:',
  119. available_1d)
  120. if dim == 2 and (method not in available_2d):
  121. raise ValueError('Invalid method, Available options are:',
  122. available_2d)
  123. if dim == 3 and (method not in available_3d):
  124. raise ValueError('Invalid method, Available options are:',
  125. available_3d)
  126. if method in kernel_methods:
  127. input_array = np.zeros((len(lfp), lfp[0].magnitude.shape[0]))
  128. for ii, jj in enumerate(lfp):
  129. input_array[ii, :] = jj.rescale(pq.mV).magnitude
  130. kernel_method = getattr(KCSD, method) # fetch the class 'KCSD1D'
  131. lambdas = kwargs.pop('lambdas', None)
  132. Rs = kwargs.pop('Rs', None)
  133. k = kernel_method(np.array(coordinates), input_array.T, **kwargs)
  134. if process_estimate:
  135. k.cross_validate(lambdas, Rs)
  136. estm_csd = k.values()
  137. estm_csd = np.rollaxis(estm_csd, -1, 0)
  138. output = neo.AnalogSignal(estm_csd * pq.uA / pq.mm**3,
  139. t_start=lfp.t_start,
  140. sampling_rate=lfp.sampling_rate)
  141. if dim == 1:
  142. output.annotate(x_coords=k.estm_x)
  143. elif dim == 2:
  144. output.annotate(x_coords=k.estm_x, y_coords=k.estm_y)
  145. elif dim == 3:
  146. output.annotate(x_coords=k.estm_x, y_coords=k.estm_y,
  147. z_coords=k.estm_z)
  148. elif method in py_iCSD_toolbox:
  149. coordinates = np.array(coordinates) * coordinates[0].units
  150. if method in icsd_methods:
  151. try:
  152. coordinates = coordinates.rescale(kwargs['diam'].units)
  153. except KeyError: # Then why specify as a default in icsd?
  154. # All iCSD methods explicitly assume a source
  155. # diameter in contrast to the stdCSD that
  156. # implicitly assume infinite source radius
  157. raise ValueError("Parameter diam must be specified for iCSD \
  158. methods: {}".format(", ".join(icsd_methods)))
  159. if 'f_type' in kwargs:
  160. if (kwargs['f_type'] != 'identity') and \
  161. (kwargs['f_order'] is None):
  162. raise ValueError("The order of {} filter must be \
  163. specified".format(kwargs['f_type']))
  164. lfp = neo.AnalogSignal(np.asarray(lfp).T, units=lfp.units,
  165. sampling_rate=lfp.sampling_rate)
  166. csd_method = getattr(icsd, method) # fetch class from icsd.py file
  167. csd_estimator = csd_method(lfp=lfp.magnitude * lfp.units,
  168. coord_electrode=coordinates.flatten(),
  169. **kwargs)
  170. csd_pqarr = csd_estimator.get_csd()
  171. if process_estimate:
  172. csd_pqarr_filtered = csd_estimator.filter_csd(csd_pqarr)
  173. output = neo.AnalogSignal(csd_pqarr_filtered.T,
  174. t_start=lfp.t_start,
  175. sampling_rate=lfp.sampling_rate)
  176. else:
  177. output = neo.AnalogSignal(csd_pqarr.T, t_start=lfp.t_start,
  178. sampling_rate=lfp.sampling_rate)
  179. output.annotate(x_coords=coordinates)
  180. return output
  181. @deprecated_alias(ele_xx='x_positions', ele_yy='y_positions',
  182. ele_zz='z_positions', xlims='x_limits', ylims='y_limits',
  183. zlims='z_limits', res='resolution')
  184. def generate_lfp(csd_profile, x_positions, y_positions=None, z_positions=None,
  185. x_limits=[0., 1.], y_limits=[0., 1.], z_limits=[0., 1.],
  186. resolution=50):
  187. """
  188. Forward modelling for getting the potentials for testing Current Source
  189. Density (CSD).
  190. Parameters
  191. ----------
  192. csd_profile : callable
  193. A function that computes true CSD profile.
  194. Available options are (see ./csd/utility_functions.py)
  195. 1D : gauss_1d_dipole
  196. 2D : large_source_2D and small_source_2D
  197. 3D : gauss_3d_dipole
  198. x_positions : np.ndarray
  199. Positions of the x coordinates of the electrodes
  200. y_positions : np.ndarray, optional
  201. Positions of the y coordinates of the electrodes
  202. Defaults to None, use in 2D or 3D cases only
  203. z_positions : np.ndarray, optional
  204. Positions of the z coordinates of the electrodes
  205. Defaults to None, use in 3D case only
  206. x_limits : list, optional
  207. A list of [start, end].
  208. The starting spatial coordinate and the ending for integration
  209. Defaults to [0.,1.]
  210. y_limits : list, optional
  211. A list of [start, end].
  212. The starting spatial coordinate and the ending for integration
  213. Defaults to [0.,1.], use only in 2D and 3D case
  214. z_limits : list, optional
  215. A list of [start, end].
  216. The starting spatial coordinate and the ending for integration
  217. Defaults to [0.,1.], use only in 3D case
  218. resolution : int, optional
  219. The resolution of the integration
  220. Defaults to 50
  221. Returns
  222. -------
  223. LFP : neo.AnalogSignal
  224. The potentials created by the csd profile at the electrode positions.
  225. The electrode positions are attached as RecordingChannel's coordinate.
  226. """
  227. def integrate_1D(x0, csd_x, csd, h):
  228. m = np.sqrt((csd_x - x0) ** 2 + h ** 2) - abs(csd_x - x0)
  229. y = csd * m
  230. I = simps(y, csd_x)
  231. return I
  232. def integrate_2D(x, y, xlin, ylin, csd, h, X, Y):
  233. x = np.reshape(x, (1, 1, len(x)))
  234. y = np.reshape(y, (1, 1, len(y)))
  235. X = np.expand_dims(X, axis=2)
  236. Y = np.expand_dims(Y, axis=2)
  237. csd = np.expand_dims(csd, axis=2)
  238. m = np.sqrt((x - X) ** 2 + (y - Y) ** 2)
  239. np.clip(m, a_min=0.0000001, a_max=None, out=m)
  240. y = np.arcsinh(2 * h / m) * csd
  241. I = simps(y.T, ylin)
  242. F = simps(I, xlin)
  243. return F
  244. def integrate_3D(x, y, z, csd, xlin, ylin, zlin, X, Y, Z):
  245. m = np.sqrt((x - X) ** 2 + (y - Y) ** 2 + (z - Z) ** 2)
  246. np.clip(m, a_min=0.0000001, a_max=None, out=m)
  247. z = csd / m
  248. Iy = simps(np.transpose(z, (1, 0, 2)), zlin)
  249. Iy = simps(Iy, ylin)
  250. F = simps(Iy, xlin)
  251. return F
  252. dim = 1
  253. if z_positions is not None:
  254. dim = 3
  255. elif y_positions is not None:
  256. dim = 2
  257. x = np.linspace(x_limits[0], x_limits[1], resolution)
  258. sigma = 1.0
  259. h = 50.
  260. if dim == 1:
  261. chrg_x = x
  262. csd = csd_profile(chrg_x)
  263. pots = integrate_1D(x_positions, chrg_x, csd, h)
  264. pots /= 2. * sigma # eq.: 26 from Potworowski et al
  265. ele_pos = x_positions
  266. elif dim == 2:
  267. y = np.linspace(y_limits[0], y_limits[1], resolution)
  268. chrg_x = np.expand_dims(x, axis=1)
  269. chrg_y = np.expand_dims(y, axis=0)
  270. csd = csd_profile(chrg_x, chrg_y)
  271. pots = integrate_2D(x_positions, y_positions,
  272. x, y,
  273. csd, h,
  274. chrg_x, chrg_y)
  275. pots /= 2 * np.pi * sigma
  276. ele_pos = np.vstack((x_positions, y_positions)).T
  277. elif dim == 3:
  278. y = np.linspace(y_limits[0], y_limits[1], resolution)
  279. z = np.linspace(z_limits[0], z_limits[1], resolution)
  280. chrg_x, chrg_y, chrg_z = np.mgrid[
  281. x_limits[0]: x_limits[1]: np.complex(0, resolution),
  282. y_limits[0]: y_limits[1]: np.complex(0, resolution),
  283. z_limits[0]: z_limits[1]: np.complex(0, resolution)
  284. ]
  285. csd = csd_profile(chrg_x, chrg_y, chrg_z)
  286. pots = np.zeros(len(x_positions))
  287. for ii in range(len(x_positions)):
  288. pots[ii] = integrate_3D(x_positions[ii], y_positions[ii],
  289. z_positions[ii],
  290. csd,
  291. x, y, z,
  292. chrg_x, chrg_y, chrg_z)
  293. pots /= 4 * np.pi * sigma
  294. ele_pos = np.vstack((x_positions, y_positions, z_positions)).T
  295. ele_pos = ele_pos * pq.mm
  296. ch = neo.ChannelIndex(index=range(len(pots)))
  297. asig = neo.AnalogSignal(np.expand_dims(pots, axis=0),
  298. sampling_rate=pq.kHz, units='mV')
  299. ch.coordinates = ele_pos
  300. ch.analogsignals.append(asig)
  301. ch.create_relationship()
  302. return asig