Skip to content

Entrainment Analysis

The entrainment module provides tools for analyzing the entrainment of spikes and LFP signals, including phase-locking value (PLV) and pairwise phase consistency (PPC) calculations.

bmtool.analysis.entrainment.align_spike_times_with_lfp(lfp, timestamps)

the lfp xarray should have a time axis. use that to align the spike times since the lfp can start at a non-zero time after sliced. Both need to be on same fs for this to be correct.

Parameters:

Name Type Description Default
lfp DataArray

LFP data with time coordinates

required
timestamps ndarray

Array of spike timestamps

required

Returns:

Type Description
ndarray

Copy of timestamps with adjusted timestamps to align with lfp.

Source code in bmtool/analysis/entrainment.py
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
def align_spike_times_with_lfp(lfp: xr.DataArray, timestamps: np.ndarray) -> np.ndarray:
    """the lfp xarray should have a time axis. use that to align the spike times since the lfp can start at a
    non-zero time after sliced. Both need to be on same fs for this to be correct.

    Parameters
    ----------
    lfp : xarray.DataArray
        LFP data with time coordinates
    timestamps : np.ndarray
        Array of spike timestamps

    Returns
    -------
    np.ndarray
        Copy of timestamps with adjusted timestamps to align with lfp.
    """
    # print("Pairing LFP and Spike Times")
    # print(lfp.time.values)
    # print(f"LFP starts at {lfp.time.values[0]}ms")
    # need to make sure lfp and spikes have the same time axis
    # align spikes with lfp
    timestamps = timestamps[
        (timestamps >= lfp.time.values[0]) & (timestamps <= lfp.time.values[-1])
    ].copy()
    # set the time axis of the spikes to match the lfp
    # timestamps = timestamps - lfp.time.values[0]
    return timestamps

bmtool.analysis.entrainment.calculate_signal_signal_plv(signal1, signal2, fs, freq_of_interest=None, filter_method='wavelet', lowcut=None, highcut=None, bandwidth=2.0)

Calculate Phase Locking Value (PLV) between two signals using wavelet or Hilbert method.

Parameters:

Name Type Description Default
signal1 ndarray

First input signal (1D array)

required
signal2 ndarray

Second input signal (1D array, same length as signal1)

required
fs float

Sampling frequency in Hz

required
freq_of_interest float

Desired frequency for wavelet PLV calculation, required if filter_method='wavelet'

None
filter_method str

Method to use for filtering, either 'wavelet' or 'butter' (default: 'wavelet')

'wavelet'
lowcut float

Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
highcut float

Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
bandwidth float

Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)

2.0

Returns:

Type Description
ndarray

Phase Locking Value (1D array)

Source code in bmtool/analysis/entrainment.py
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
def calculate_signal_signal_plv(
    signal1: np.ndarray,
    signal2: np.ndarray,
    fs: float,
    freq_of_interest: float = None,
    filter_method: str = "wavelet",
    lowcut: float = None,
    highcut: float = None,
    bandwidth: float = 2.0,
) -> np.ndarray:
    """
    Calculate Phase Locking Value (PLV) between two signals using wavelet or Hilbert method.

    Parameters
    ----------
    signal1 : np.ndarray
        First input signal (1D array)
    signal2 : np.ndarray
        Second input signal (1D array, same length as signal1)
    fs : float
        Sampling frequency in Hz
    freq_of_interest : float, optional
        Desired frequency for wavelet PLV calculation, required if filter_method='wavelet'
    filter_method : str, optional
        Method to use for filtering, either 'wavelet' or 'butter' (default: 'wavelet')
    lowcut : float, optional
        Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    highcut : float, optional
        Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    bandwidth : float, optional
        Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)

    Returns
    -------
    np.ndarray
        Phase Locking Value (1D array)
    """
    if len(signal1) != len(signal2):
        raise ValueError("Input signals must have the same length.")

    if filter_method == "wavelet":
        if freq_of_interest is None:
            raise ValueError("freq_of_interest must be provided for the wavelet method.")

        # Apply CWT to both signals
        theta1 = wavelet_filter(x=signal1, freq=freq_of_interest, fs=fs, bandwidth=bandwidth)
        theta2 = wavelet_filter(x=signal2, freq=freq_of_interest, fs=fs, bandwidth=bandwidth)

    elif filter_method == "butter":
        if lowcut is None or highcut is None:
            print(
                "Lowcut and/or highcut were not defined, signal will not be filtered and will just take Hilbert transform for PLV calculation"
            )

        if lowcut and highcut:
            # Bandpass filter and get the analytic signal using the Hilbert transform
            filtered_signal1 = butter_bandpass_filter(
                data=signal1, lowcut=lowcut, highcut=highcut, fs=fs
            )
            filtered_signal2 = butter_bandpass_filter(
                data=signal2, lowcut=lowcut, highcut=highcut, fs=fs
            )
            # Get phase using the Hilbert transform
            theta1 = signal.hilbert(filtered_signal1)
            theta2 = signal.hilbert(filtered_signal2)
        else:
            # Get phase using the Hilbert transform without filtering
            theta1 = signal.hilbert(signal1)
            theta2 = signal.hilbert(signal2)

    else:
        raise ValueError("Invalid method. Choose 'wavelet' or 'butter'.")

    # Calculate phase difference
    phase_diff = np.angle(theta1) - np.angle(theta2)

    # Calculate PLV from standard equation from Measuring phase synchrony in brain signals(1999)
    plv = np.abs(np.mean(np.exp(1j * phase_diff), axis=-1))

    return plv

bmtool.analysis.entrainment.calculate_spike_lfp_plv(spike_times=None, lfp_data=None, spike_fs=None, lfp_fs=None, filter_method='butter', freq_of_interest=None, lowcut=None, highcut=None, bandwidth=2.0, filtered_lfp_phase=None)

Calculate spike-lfp unbiased phase locking value

Parameters:

Name Type Description Default
spike_times ndarray

Array of spike times

None
lfp_data ndarray

Local field potential time series data. Not required if filtered_lfp_phase is provided.

None
spike_fs float

Sampling frequency in Hz of the spike times, only needed if spike times and LFP have different sampling rates

None
lfp_fs float

Sampling frequency in Hz of the LFP data

None
filter_method str

Method to use for filtering, either 'wavelet' or 'butter' (default: 'butter')

'butter'
freq_of_interest float

Desired frequency for wavelet phase extraction, required if filter_method='wavelet'

None
lowcut float

Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
highcut float

Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
bandwidth float

Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)

2.0
filtered_lfp_phase ndarray

Pre-computed instantaneous phase of the filtered LFP. If provided, the function will skip the filtering step.

None

Returns:

Type Description
float

Phase Locking Value (unbiased)

Source code in bmtool/analysis/entrainment.py
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
def calculate_spike_lfp_plv(
    spike_times: np.ndarray = None,
    lfp_data=None,
    spike_fs: float = None,
    lfp_fs: float = None,
    filter_method: str = "butter",
    freq_of_interest: float = None,
    lowcut: float = None,
    highcut: float = None,
    bandwidth: float = 2.0,
    filtered_lfp_phase: np.ndarray = None,
) -> float:
    """
    Calculate spike-lfp unbiased phase locking value

    Parameters
    ----------
    spike_times : np.ndarray
        Array of spike times
    lfp_data : np.ndarray
        Local field potential time series data. Not required if filtered_lfp_phase is provided.
    spike_fs : float, optional
        Sampling frequency in Hz of the spike times, only needed if spike times and LFP have different sampling rates
    lfp_fs : float
        Sampling frequency in Hz of the LFP data
    filter_method : str, optional
        Method to use for filtering, either 'wavelet' or 'butter' (default: 'butter')
    freq_of_interest : float, optional
        Desired frequency for wavelet phase extraction, required if filter_method='wavelet'
    lowcut : float, optional
        Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    highcut : float, optional
        Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    bandwidth : float, optional
        Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)
    filtered_lfp_phase : np.ndarray, optional
        Pre-computed instantaneous phase of the filtered LFP. If provided, the function will skip the filtering step.

    Returns
    -------
    float
        Phase Locking Value (unbiased)
    """

    spike_phases = _get_spike_phases(
        spike_times=spike_times,
        lfp_data=lfp_data,
        spike_fs=spike_fs,
        lfp_fs=lfp_fs,
        filter_method=filter_method,
        freq_of_interest=freq_of_interest,
        lowcut=lowcut,
        highcut=highcut,
        bandwidth=bandwidth,
        filtered_lfp_phase=filtered_lfp_phase,
    )

    if len(spike_phases) <= 1:
        return 0

    # Number of spikes
    N = len(spike_phases)

    # Convert phases to unit vectors in the complex plane
    unit_vectors = np.exp(1j * spike_phases)

    # Sum of all unit vectors (resultant vector)
    resultant_vector = np.sum(unit_vectors)

    # Calculate plv^2 * N
    plv2n = (resultant_vector * resultant_vector.conjugate()).real / N  # plv^2 * N
    plv = (plv2n / N) ** 0.5
    ppc = (plv2n - 1) / (N - 1)  # ppc = (plv^2 * N - 1) / (N - 1)
    plv_unbiased = np.fmax(ppc, 0.0) ** 0.5  # ensure non-negative

    return plv_unbiased

bmtool.analysis.entrainment.calculate_ppc(spike_times=None, lfp_data=None, spike_fs=None, lfp_fs=None, filter_method='wavelet', freq_of_interest=None, lowcut=None, highcut=None, bandwidth=2.0, ppc_method='numpy', filtered_lfp_phase=None)

Calculate Pairwise Phase Consistency (PPC) between spike times and LFP signal. Based on https://www.sciencedirect.com/science/article/pii/S1053811910000959

Parameters:

Name Type Description Default
spike_times ndarray

Array of spike times

None
lfp_data ndarray

Local field potential time series data. Not required if filtered_lfp_phase is provided.

None
spike_fs float

Sampling frequency in Hz of the spike times, only needed if spike times and LFP have different sampling rates

None
lfp_fs float

Sampling frequency in Hz of the LFP data

None
filter_method str

Method to use for filtering, either 'wavelet' or 'butter' (default: 'wavelet')

'wavelet'
freq_of_interest float

Desired frequency for wavelet phase extraction, required if filter_method='wavelet'

None
lowcut float

Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
highcut float

Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
bandwidth float

Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)

2.0
ppc_method str

Algorithm to use for PPC calculation: 'numpy', 'numba', or 'gpu' (default: 'numpy')

'numpy'
filtered_lfp_phase ndarray

Pre-computed instantaneous phase of the filtered LFP. If provided, the function will skip the filtering step.

None

Returns:

Type Description
float

Pairwise Phase Consistency value

Source code in bmtool/analysis/entrainment.py
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
def calculate_ppc(
    spike_times: np.ndarray = None,
    lfp_data=None,
    spike_fs: float = None,
    lfp_fs: float = None,
    filter_method: str = "wavelet",
    freq_of_interest: float = None,
    lowcut: float = None,
    highcut: float = None,
    bandwidth: float = 2.0,
    ppc_method: str = "numpy",
    filtered_lfp_phase: np.ndarray = None,
) -> float:
    """
    Calculate Pairwise Phase Consistency (PPC) between spike times and LFP signal.
    Based on https://www.sciencedirect.com/science/article/pii/S1053811910000959

    Parameters
    ----------
    spike_times : np.ndarray
        Array of spike times
    lfp_data : np.ndarray
        Local field potential time series data. Not required if filtered_lfp_phase is provided.
    spike_fs : float, optional
        Sampling frequency in Hz of the spike times, only needed if spike times and LFP have different sampling rates
    lfp_fs : float
        Sampling frequency in Hz of the LFP data
    filter_method : str, optional
        Method to use for filtering, either 'wavelet' or 'butter' (default: 'wavelet')
    freq_of_interest : float, optional
        Desired frequency for wavelet phase extraction, required if filter_method='wavelet'
    lowcut : float, optional
        Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    highcut : float, optional
        Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    bandwidth : float, optional
        Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)
    ppc_method : str, optional
        Algorithm to use for PPC calculation: 'numpy', 'numba', or 'gpu' (default: 'numpy')
    filtered_lfp_phase : np.ndarray, optional
        Pre-computed instantaneous phase of the filtered LFP. If provided, the function will skip the filtering step.

    Returns
    -------
    float
        Pairwise Phase Consistency value
    """

    spike_phases = _get_spike_phases(
        spike_times=spike_times,
        lfp_data=lfp_data,
        spike_fs=spike_fs,
        lfp_fs=lfp_fs,
        filter_method=filter_method,
        freq_of_interest=freq_of_interest,
        lowcut=lowcut,
        highcut=highcut,
        bandwidth=bandwidth,
        filtered_lfp_phase=filtered_lfp_phase,
    )

    if len(spike_phases) <= 1:
        return 0

    n_spikes = len(spike_phases)

    # Calculate PPC (Pairwise Phase Consistency)
    # Explicit calculation of pairwise phase consistency
    # Vectorized computation for efficiency
    if ppc_method == "numpy":
        i, j = np.triu_indices(n_spikes, k=1)
        phase_diff = spike_phases[i] - spike_phases[j]
        sum_cos_diff = np.sum(np.cos(phase_diff))
        ppc = (2 / (n_spikes * (n_spikes - 1))) * sum_cos_diff
    elif ppc_method == "numba":
        ppc = _ppc_parallel_numba(spike_phases)
    elif ppc_method == "gpu":
        ppc = _ppc_gpu(spike_phases)
    else:
        raise ValueError("Please use a supported ppc method currently that is numpy, numba or gpu")
    return ppc

bmtool.analysis.entrainment.calculate_ppc2(spike_times=None, lfp_data=None, spike_fs=None, lfp_fs=None, filter_method='wavelet', freq_of_interest=None, lowcut=None, highcut=None, bandwidth=2.0, filtered_lfp_phase=None)

-----------------------------------------------------------------------------

PPC2 Calculation (Vinck et al., 2010)

-----------------------------------------------------------------------------

Equation(Original):

PPC = (2 / (n * (n - 1))) * sum(cos(φ_i - φ_j) for all i < j)

Optimized Formula (Algebraically Equivalent):

PPC = (|sum(e^(i*φ_j))|^2 - n) / (n * (n - 1))

-----------------------------------------------------------------------------

Parameters:

Name Type Description Default
spike_times ndarray

Array of spike times

None
lfp_data ndarray

Local field potential time series data. Not required if filtered_lfp_phase is provided.

None
spike_fs float

Sampling frequency in Hz of the spike times, only needed if spike times and LFP have different sampling rates

None
lfp_fs float

Sampling frequency in Hz of the LFP data

None
filter_method str

Method to use for filtering, either 'wavelet' or 'butter' (default: 'wavelet')

'wavelet'
freq_of_interest float

Desired frequency for wavelet phase extraction, required if filter_method='wavelet'

None
lowcut float

Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
highcut float

Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
bandwidth float

Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)

2.0
filtered_lfp_phase ndarray

Pre-computed instantaneous phase of the filtered LFP. If provided, the function will skip the filtering step.

None

Returns:

Type Description
float

Pairwise Phase Consistency 2 (PPC2) value

Source code in bmtool/analysis/entrainment.py
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
def calculate_ppc2(
    spike_times: np.ndarray = None,
    lfp_data=None,
    spike_fs: float = None,
    lfp_fs: float = None,
    filter_method: str = "wavelet",
    freq_of_interest: float = None,
    lowcut: float = None,
    highcut: float = None,
    bandwidth: float = 2.0,
    filtered_lfp_phase: np.ndarray = None,
) -> float:
    """
    # -----------------------------------------------------------------------------
    # PPC2 Calculation (Vinck et al., 2010)
    # -----------------------------------------------------------------------------
    # Equation(Original):
    #   PPC = (2 / (n * (n - 1))) * sum(cos(φ_i - φ_j) for all i < j)
    # Optimized Formula (Algebraically Equivalent):
    #   PPC = (|sum(e^(i*φ_j))|^2 - n) / (n * (n - 1))
    # -----------------------------------------------------------------------------

    Parameters
    ----------
    spike_times : np.ndarray
        Array of spike times
    lfp_data : np.ndarray
        Local field potential time series data. Not required if filtered_lfp_phase is provided.
    spike_fs : float, optional
        Sampling frequency in Hz of the spike times, only needed if spike times and LFP have different sampling rates
    lfp_fs : float
        Sampling frequency in Hz of the LFP data
    filter_method : str, optional
        Method to use for filtering, either 'wavelet' or 'butter' (default: 'wavelet')
    freq_of_interest : float, optional
        Desired frequency for wavelet phase extraction, required if filter_method='wavelet'
    lowcut : float, optional
        Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    highcut : float, optional
        Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    bandwidth : float, optional
        Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)
    filtered_lfp_phase : np.ndarray, optional
        Pre-computed instantaneous phase of the filtered LFP. If provided, the function will skip the filtering step.

    Returns
    -------
    float
        Pairwise Phase Consistency 2 (PPC2) value
    """

    spike_phases = _get_spike_phases(
        spike_times=spike_times,
        lfp_data=lfp_data,
        spike_fs=spike_fs,
        lfp_fs=lfp_fs,
        filter_method=filter_method,
        freq_of_interest=freq_of_interest,
        lowcut=lowcut,
        highcut=highcut,
        bandwidth=bandwidth,
        filtered_lfp_phase=filtered_lfp_phase,
    )

    if len(spike_phases) <= 1:
        return 0

    # Calculate PPC2 according to Vinck et al. (2010), Equation 6
    n = len(spike_phases)

    # Convert phases to unit vectors in the complex plane
    unit_vectors = np.exp(1j * spike_phases)

    # Calculate the resultant vector
    resultant_vector = np.sum(unit_vectors)

    # PPC2 = (|∑(e^(i*φ_j))|² - n) / (n * (n - 1))
    ppc2 = (np.abs(resultant_vector) ** 2 - n) / (n * (n - 1))

    return ppc2

bmtool.analysis.entrainment.calculate_entrainment_per_cell(spike_df=None, lfp_data=None, filter_method='wavelet', pop_names=None, entrainment_method='plv', lowcut=None, highcut=None, spike_fs=None, lfp_fs=None, bandwidth=2, freqs=None, ppc_method='numpy')

Calculate neural entrainment (PPC, PLV) per neuron (cell) for specified frequencies across different populations.

This function computes the entrainment metrics for each neuron within the specified populations based on their spike times and the provided LFP signal. It returns a nested dictionary structure containing the entrainment values organized by population, node ID, and frequency.

Parameters:

Name Type Description Default
spike_df DataFrame

DataFrame containing spike data with columns 'pop_name', 'node_ids', and 'timestamps'

None
lfp_data ndarray

Local field potential (LFP) time series data

None
filter_method str

Method to use for filtering, either 'wavelet' or 'butter' (default: 'wavelet')

'wavelet'
entrainment_method str

Method to use for entrainment calculation, either 'plv', 'ppc', or 'ppc2' (default: 'plv')

'plv'
lowcut float

Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
highcut float

Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
spike_fs float

Sampling frequency of the spike times in Hz

None
lfp_fs float

Sampling frequency of the LFP signal in Hz

None
bandwidth float

Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)

2
ppc_method str

Algorithm to use for PPC calculation: 'numpy', 'numba', or 'gpu' (default: 'numpy')

'numpy'
pop_names List[str]

List of population names to analyze

None
freqs List[float]

List of frequencies (in Hz) at which to calculate entrainment

None

Returns:

Type Description
Dict[str, Dict[int, Dict[float, float]]]

Nested dictionary where the structure is: { population_name: { node_id: { frequency: entrainment value } } } Entrainment values are floats representing the metric (PPC, PLV) at each frequency

Source code in bmtool/analysis/entrainment.py
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
def calculate_entrainment_per_cell(
    spike_df: pd.DataFrame = None,
    lfp_data: np.ndarray = None,
    filter_method: str = "wavelet",
    pop_names: List[str] = None,
    entrainment_method: str = "plv",
    lowcut: float = None,
    highcut: float = None,
    spike_fs: float = None,
    lfp_fs: float = None,
    bandwidth: float = 2,
    freqs: List[float] = None,
    ppc_method: str = "numpy",
) -> Dict[str, Dict[int, Dict[float, float]]]:
    """
    Calculate neural entrainment (PPC, PLV) per neuron (cell) for specified frequencies across different populations.

    This function computes the entrainment metrics for each neuron within the specified populations based on their spike times
    and the provided LFP signal. It returns a nested dictionary structure containing the entrainment values
    organized by population, node ID, and frequency.

    Parameters
    ----------
    spike_df : pd.DataFrame
        DataFrame containing spike data with columns 'pop_name', 'node_ids', and 'timestamps'
    lfp_data : np.ndarray
        Local field potential (LFP) time series data
    filter_method : str, optional
        Method to use for filtering, either 'wavelet' or 'butter' (default: 'wavelet')
    entrainment_method : str, optional
        Method to use for entrainment calculation, either 'plv', 'ppc', or 'ppc2' (default: 'plv')
    lowcut : float, optional
        Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    highcut : float, optional
        Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    spike_fs : float
        Sampling frequency of the spike times in Hz
    lfp_fs : float
        Sampling frequency of the LFP signal in Hz
    bandwidth : float, optional
        Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)
    ppc_method : str, optional
        Algorithm to use for PPC calculation: 'numpy', 'numba', or 'gpu' (default: 'numpy')
    pop_names : List[str]
        List of population names to analyze
    freqs : List[float]
        List of frequencies (in Hz) at which to calculate entrainment

    Returns
    -------
    Dict[str, Dict[int, Dict[float, float]]]
        Nested dictionary where the structure is:
        {
            population_name: {
                node_id: {
                    frequency: entrainment value
                }
            }
        }
        Entrainment values are floats representing the metric (PPC, PLV) at each frequency
    """
    # pre filter lfp to speed up calculate of entrainment
    filtered_lfp_phases = {}
    for freq in range(len(freqs)):
        phase = get_lfp_phase(
            lfp_data=lfp_data,
            freq_of_interest=freqs[freq],
            fs=lfp_fs,
            filter_method=filter_method,
            lowcut=lowcut,
            highcut=highcut,
            bandwidth=bandwidth,
        )
        filtered_lfp_phases[freqs[freq]] = phase

    entrainment_dict = {}
    for pop in pop_names:
        skip_count = 0
        pop_spikes = spike_df[spike_df["pop_name"] == pop]
        nodes = sorted(pop_spikes["node_ids"].unique())  # sort so all nodes are processed in order
        entrainment_dict[pop] = {}
        print(f"Processing {pop} population")
        for node in tqdm(nodes):
            node_spikes = pop_spikes[pop_spikes["node_ids"] == node]

            # Skip nodes with less than or equal to 1 spike
            if len(node_spikes) <= 1:
                skip_count += 1
                continue

            entrainment_dict[pop][node] = {}
            for freq in freqs:
                # Calculate entrainment based on the selected method using the pre-filtered phases
                if entrainment_method == "plv":
                    entrainment_dict[pop][node][freq] = calculate_spike_lfp_plv(
                        node_spikes["timestamps"].values,
                        lfp_data,
                        spike_fs=spike_fs,
                        lfp_fs=lfp_fs,
                        freq_of_interest=freq,
                        bandwidth=bandwidth,
                        lowcut=lowcut,
                        highcut=highcut,
                        filter_method=filter_method,
                        filtered_lfp_phase=filtered_lfp_phases[freq],
                    )
                elif entrainment_method == "ppc2":
                    entrainment_dict[pop][node][freq] = calculate_ppc2(
                        node_spikes["timestamps"].values,
                        lfp_data,
                        spike_fs=spike_fs,
                        lfp_fs=lfp_fs,
                        freq_of_interest=freq,
                        bandwidth=bandwidth,
                        lowcut=lowcut,
                        highcut=highcut,
                        filter_method=filter_method,
                        filtered_lfp_phase=filtered_lfp_phases[freq],
                    )
                elif entrainment_method == "ppc":
                    entrainment_dict[pop][node][freq] = calculate_ppc(
                        node_spikes["timestamps"].values,
                        lfp_data,
                        spike_fs=spike_fs,
                        lfp_fs=lfp_fs,
                        freq_of_interest=freq,
                        bandwidth=bandwidth,
                        lowcut=lowcut,
                        highcut=highcut,
                        filter_method=filter_method,
                        ppc_method=ppc_method,
                        filtered_lfp_phase=filtered_lfp_phases[freq],
                    )

        print(
            f"Calculated {entrainment_method.upper()} for {pop} population with {len(nodes)-skip_count} valid cells, skipped {skip_count} cells for lack of spikes"
        )

    return entrainment_dict

bmtool.analysis.entrainment.calculate_spike_rate_power_correlation(spike_rate, lfp_data, fs, pop_names, filter_method='wavelet', bandwidth=2.0, lowcut=None, highcut=None, freq_range=(10, 100), freq_step=5, type_name='raw')

Calculate correlation between population spike rates (xarray) and LFP power across frequencies.

Parameters:

Name Type Description Default
spike_rate DataArray

Population spike rates with dimensions (time, population[, type])

required
lfp_data ndarray

LFP data

required
fs float

Sampling frequency

required
pop_names list

List of population names to analyze

required
filter_method str

Filtering method to use, either 'wavelet' or 'butter' (default: 'wavelet')

'wavelet'
bandwidth float

Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)

2.0
lowcut float

Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
highcut float

Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'

None
freq_range tuple

Min and max frequency to analyze (default: (10, 100))

(10, 100)
freq_step float

Step size for frequency analysis (default: 5)

5
type_name str

Which type of spike rate to use if 'type' dimension exists (default: 'raw')

'raw'

Returns:

Name Type Description
correlation_results dict

Dictionary with correlation results for each population and frequency

frequencies array

Array of frequencies analyzed

Source code in bmtool/analysis/entrainment.py
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
def calculate_spike_rate_power_correlation(
    spike_rate: xr.DataArray,
    lfp_data: np.ndarray,
    fs: float,
    pop_names: list,
    filter_method: str = "wavelet",
    bandwidth: float = 2.0,
    lowcut: float = None,
    highcut: float = None,
    freq_range: tuple = (10, 100),
    freq_step: float = 5,
    type_name: str = "raw",  # 'raw' or 'smoothed'
):
    """
    Calculate correlation between population spike rates (xarray) and LFP power across frequencies.

    Parameters
    ----------
    spike_rate : xr.DataArray
        Population spike rates with dimensions (time, population[, type])
    lfp_data : np.ndarray
        LFP data
    fs : float
        Sampling frequency
    pop_names : list
        List of population names to analyze
    filter_method : str, optional
        Filtering method to use, either 'wavelet' or 'butter' (default: 'wavelet')
    bandwidth : float, optional
        Bandwidth parameter for wavelet filter when method='wavelet' (default: 2.0)
    lowcut : float, optional
        Lower frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    highcut : float, optional
        Upper frequency bound (Hz) for butterworth bandpass filter, required if filter_method='butter'
    freq_range : tuple, optional
        Min and max frequency to analyze (default: (10, 100))
    freq_step : float, optional
        Step size for frequency analysis (default: 5)
    type_name : str, optional
        Which type of spike rate to use if 'type' dimension exists (default: 'raw')

    Returns
    -------
    correlation_results : dict
        Dictionary with correlation results for each population and frequency
    frequencies : array
        Array of frequencies analyzed
    """
    frequencies = np.arange(freq_range[0], freq_range[1] + 1, freq_step)
    correlation_results = {pop: {} for pop in pop_names}

    # Calculate power at each frequency band using specified filter
    power_by_freq = {}
    for freq in frequencies:
        power_by_freq[freq] = get_lfp_power(
            lfp_data, freq, fs, filter_method, lowcut=lowcut, highcut=highcut, bandwidth=bandwidth
        )

    # For each population, extract the correct spike rate
    for pop in pop_names:
        # If 'type' dimension exists, select the type
        if "type" in spike_rate.dims:
            pop_rate = spike_rate.sel(population=pop, type=type_name).values
        else:
            pop_rate = spike_rate.sel(population=pop).values

        # Calculate correlation with power at each frequency
        for freq in frequencies:
            lfp_power = power_by_freq[freq]
            # Ensure lengths match
            min_len = min(len(pop_rate), len(lfp_power))
            if len(pop_rate) != len(lfp_power):
                print(f"Warning: Length mismatch for {pop} at {freq} Hz, truncating to {min_len}")
            corr, p_val = stats.spearmanr(pop_rate[:min_len], lfp_power[:min_len])
            correlation_results[pop][freq] = {"correlation": corr, "p_value": p_val}

    return correlation_results, frequencies

bmtool.analysis.entrainment.get_spikes_in_cycle(spike_df, lfp_data, spike_fs=1000, lfp_fs=400, filter_method='butter', lowcut=None, highcut=None, bandwidth=2.0, freq_of_interest=None)

Analyze spike timing relative to oscillation phases.

Parameters:

spike_df : pd.DataFrame lfp_data : np.array Raw LFP signal fs : float Sampling frequency of LFP in Hz gamma_band : tuple Lower and upper bounds of gamma frequency band in Hz

Returns:

phase_data : dict Dictionary containing phase values for each spike and neuron population

Source code in bmtool/analysis/entrainment.py
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
def get_spikes_in_cycle(
    spike_df,
    lfp_data,
    spike_fs=1000,
    lfp_fs=400,
    filter_method="butter",
    lowcut=None,
    highcut=None,
    bandwidth=2.0,
    freq_of_interest=None,
):
    """
    Analyze spike timing relative to oscillation phases.

    Parameters:
    -----------
    spike_df : pd.DataFrame
    lfp_data : np.array
        Raw LFP signal
    fs : float
        Sampling frequency of LFP in Hz
    gamma_band : tuple
        Lower and upper bounds of gamma frequency band in Hz

    Returns:
    --------
    phase_data : dict
        Dictionary containing phase values for each spike and neuron population
    """
    phase = get_lfp_phase(
        lfp_data=lfp_data,
        fs=lfp_fs,
        filter_method=filter_method,
        lowcut=lowcut,
        highcut=highcut,
        bandwidth=bandwidth,
        freq_of_interest=freq_of_interest,
    )

    # Get unique neuron populations
    neuron_pops = spike_df["pop_name"].unique()

    # Get the phase at each spike time for each neuron population
    phase_data = {}

    for pop in neuron_pops:
        # Get spike times for this population
        pop_spikes = spike_df[spike_df["pop_name"] == pop]["timestamps"].values

        # Convert spike times to sample indices
        spike_times_seconds = pop_spikes / spike_fs

        # Then convert from seconds to samples at the new sampling rate
        spike_indices = np.round(spike_times_seconds * lfp_fs).astype(int)

        # Ensure spike times are within LFP data range
        valid_indices = (spike_indices >= 0) & (spike_indices < len(phase))

        if np.any(valid_indices):
            valid_samples = spike_indices[valid_indices]
            phase_data[pop] = phase[valid_samples]

    return phase_data