netneurotools.stats.gen_spinsamples

netneurotools.stats.gen_spinsamples(coords, hemiid, n_rotate=1000, check_duplicates=True, method='original', exact=False, seed=None, verbose=False, return_cost=False)[source]

Return a resampling array for coords obtained from rotations / spins.

Using the method initially proposed in [ST1] (and later modified + updated based on findings in [ST2] and [ST3]), this function applies random rotations to the user-supplied coords in order to generate a resampling array that preserves its spatial embedding. Rotations are generated for one hemisphere and mirrored for the other (see hemiid for more information).

Due to irregular sampling of coords and the randomness of the rotations it is possible that some “rotations” may resample with replacement (i.e., will not be a true permutation). The likelihood of this can be reduced by either increasing the sampling density of coords or changing the method parameter (see Notes for more information on the latter).

Parameters:
  • coords ((N, 3) array_like) – X, Y, Z coordinates of N nodes/parcels/regions/vertices defined on a sphere

  • hemiid ((N,) array_like) – Array denoting hemisphere designation of coordinates in coords, where values should be {0, 1} denoting the different hemispheres. Rotations are generated for one hemisphere and mirrored across the y-axis for the other hemisphere.

  • n_rotate (int, optional) – Number of rotations to generate. Default: 1000

  • check_duplicates (bool, optional) – Whether to check for and attempt to avoid duplicate resamplings. A warnings will be raised if duplicates cannot be avoided. Setting to True may increase the runtime of this function! Default: True

  • method ({'original', 'vasa', 'hungarian'}, optional) – Method by which to match non- and rotated coordinates. Specifying ‘original’ will use the method described in [ST1]. Specfying ‘vasa’ will use the method described in [ST4]. Specfying ‘hungarian’ will use the Hungarian algorithm to minimize the global cost of reassignment (will dramatically increase runtime). Default: ‘original’

  • seed ({int, np.random.RandomState instance, None}, optional) – Seed for random number generation. Default: None

  • verbose (bool, optional) – Whether to print occasional status messages. Default: False

  • return_cost (bool, optional) – Whether to return cost array (specified as Euclidean distance) for each coordinate for each rotation Default: True

Returns:

  • spinsamples ((N, n_rotate) numpy.ndarray) – Resampling matrix to use in permuting data based on supplied coords.

  • cost ((N, n_rotate,) numpy.ndarray) – Cost (specified as Euclidean distance) of re-assigning each coordinate for every rotation in spinsamples. Only provided if return_cost is True.

Notes

By default, this function uses the minimum Euclidean distance between the original coordinates and the new, rotated coordinates to generate a resampling array after each spin. Unfortunately, this can (with some frequency) lead to multiple coordinates being re-assigned the same value:

>>> from netneurotools import stats as nnstats
>>> coords = [[0, 0, 1], [1, 0, 0], [0, 0, 1], [1, 0, 0]]
>>> hemi = [0, 0, 1, 1]
>>> nnstats.gen_spinsamples(coords, hemi, n_rotate=1, seed=1,
...                         method='original', check_duplicates=False)
array([[0],
       [0],
       [2],
       [3]])

While this is reasonable in most circumstances, if you feel incredibly strongly about having a perfect “permutation” (i.e., all indices appear once and exactly once in the resampling), you can set the method parameter to either ‘vasa’ or ‘hungarian’:

>>> nnstats.gen_spinsamples(coords, hemi, n_rotate=1, seed=1,
...                         method='vasa', check_duplicates=False)
array([[1],
       [0],
       [2],
       [3]])
>>> nnstats.gen_spinsamples(coords, hemi, n_rotate=1, seed=1,
...                         method='hungarian', check_duplicates=False)
array([[0],
       [1],
       [2],
       [3]])

Note that setting this parameter may increase the runtime of the function (especially for method=’hungarian’). Refer to [ST1] for information on why the default (i.e., exact set to False) suffices in most cases.

For the original MATLAB implementation of this function refer to [ST5].

References

[ST1] (1,2,3)

Alexander-Bloch, A., Shou, H., Liu, S., Satterthwaite, T. D., Glahn, D. C., Shinohara, R. T., Vandekar, S. N., & Raznahan, A. (2018). On testing for spatial correspondence between maps of human brain structure and function. NeuroImage, 178, 540-51.

[ST2]

Blaser, R., & Fryzlewicz, P. (2016). Random Rotation Ensembles. Journal of Machine Learning Research, 17(4), 1–26.

[ST3]

Lefèvre, J., Pepe, A., Muscato, J., De Guio, F., Girard, N., Auzias, G., & Germanaud, D. (2018). SPANOL (SPectral ANalysis of Lobes): A Spectral Clustering Framework for Individual and Group Parcellation of Cortical Surfaces in Lobes. Frontiers in Neuroscience, 12, 354.

[ST4]

Váša, F., Seidlitz, J., Romero-Garcia, R., Whitaker, K. J., Rosenthal, G., Vértes, P. E., … & Jones, P. B. (2018). Adolescent tuning of association cortex in human structural brain networks. Cerebral Cortex, 28(1), 281-294.