Skip to content

Signal

Source code in neurolib/utils/signal.py
 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
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
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
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
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
413
414
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
495
496
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
636
637
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
class Signal:
    name = ""
    label = ""
    signal_type = ""
    unit = ""
    description = ""
    _copy_attributes = [
        "name",
        "label",
        "signal_type",
        "unit",
        "description",
        "process_steps",
    ]
    PROCESS_STEPS_KEY = "process_steps"

    @classmethod
    def from_model_output(cls, model, group="", time_in_ms=True):
        """
        Initial Signal from modelling output.
        """
        assert isinstance(model, Model)
        return cls(model.xr(group=group), time_in_ms=time_in_ms)

    @classmethod
    def from_file(cls, filename):
        """
        Load signal from saved file.

        :param filename: filename for the Signal
        :type filename: str
        """
        if not filename.endswith(NC_EXT):
            filename += NC_EXT
        # load NC file
        xarray = xr.load_dataarray(filename)
        # init class
        signal = cls(xarray)
        # if nc file has attributes, copy them to signal class
        if xarray.attrs:
            process_steps = []
            for k, v in xarray.attrs.items():
                if cls.PROCESS_STEPS_KEY in k:
                    idx = int(k[len(cls.PROCESS_STEPS_KEY) + 1 :])
                    process_steps.insert(idx, v)
                else:
                    setattr(signal, k, v)
        else:
            logging.warning("No metadata found, setting empty...")
            process_steps = [f"raw {signal.signal_type} signal: {signal.start_time}--" f"{signal.end_time}s"]
        setattr(signal, cls.PROCESS_STEPS_KEY, process_steps)
        return signal

    def __init__(self, data, time_in_ms=False):
        """
        :param data: data for the signal, assumes time dimension with time in seconds
        :type data: xr.DataArray
        :param time_in_ms: whether time dimension is in ms
        :type time_in_ms: bool
        """
        assert isinstance(data, xr.DataArray)
        data = deepcopy(data)
        assert "time" in data.dims, "DataArray must have time axis"
        if time_in_ms:
            data["time"] = data["time"] / 1000.0
        data["time"] = np.around(data["time"], 6)
        self.data = data
        # assert time dimension is last
        self.data = self.data.transpose(*(self.dims_not_time + ["time"]))
        # compute dt and sampling frequency
        self.dt = np.around(np.diff(data.time).mean(), 6)
        self.sampling_frequency = 1.0 / self.dt
        self.process_steps = [f"raw {self.signal_type} signal: {self.start_time}--{self.end_time}s"]

    def __str__(self):
        """
        String representation.
        """
        return (
            f"{self.name} representing {self.signal_type} signal with unit of "
            f"{self.unit} with user-provided description: `{self.description}`"
            f". Shape of the signal is {self.shape} with dimensions "
            f"{self.data.dims}. Signal starts at {self.start_time} and ends at "
            f"{self.end_time}."
        )

    def __repr__(self):
        """
        Representation.
        """
        return self.__str__()

    def __eq__(self, other):
        """
        Comparison operator.

        :param other: other `Signal` to compare with
        :type other: `Signal`
        :return: whether two `Signals` are the same
        :rtype: bool
        """
        assert isinstance(other, Signal)
        # assert data are the same
        try:
            xr.testing.assert_allclose(self.data, other.data)
            eq = True
        except AssertionError:
            eq = False
        # check attributes, but if not equal, only warn the user
        for attr in self._copy_attributes:
            if getattr(self, attr) != getattr(other, attr):
                logging.warning(f"`{attr}` not equal between signals.")
        return eq

    def __getitem__(self, pos):
        """
        Get item selects in output dimension.
        """
        add_steps = [f"select `{pos}` output"]
        return self.__constructor__(self.data.sel(output=pos)).__finalize__(self, add_steps)

    def __finalize__(self, other, add_steps=None):
        """
        Copy attributes from other to self. Used when constructing class
        instance with different data, but same metadata.

        :param other: other instance of `Signal`
        :type other: `Signal`
        :param add_steps: add steps to preprocessing
        :type add_steps: list|None
        """
        assert isinstance(other, Signal)
        for attr in self._copy_attributes:
            setattr(self, attr, deepcopy(getattr(other, attr)))
        if add_steps is not None:
            self.process_steps += add_steps
        return self

    @property
    def __constructor__(self):
        """
        Return constructor, so that each child class would initiate a new
        instance of the correct class, i.e. first in the method resolution
        order.
        """
        return self.__class__.mro()[0]

    def _write_attrs_to_xr(self):
        """
        Copy attributes to xarray before saving.
        """
        # write attributes to xarray
        for attr in self._copy_attributes:
            value = getattr(self, attr)
            # if list need to unwrap
            if isinstance(value, (list, tuple)):
                for idx, val in enumerate(value):
                    self.data.attrs[f"{attr}_{idx}"] = val
            else:
                self.data.attrs[attr] = deepcopy(value)

    def save(self, filename):
        """
        Save signal.

        :param filename: filename to save, currently saves to netCDF file, which is natively supported by xarray
        :type filename: str
        """
        self._write_attrs_to_xr()
        if not filename.endswith(NC_EXT):
            filename += NC_EXT
        self.data.to_netcdf(filename)

    def iterate(self, return_as="signal"):
        """
        Return iterator over columns, so univariate measures can be computed
        per column. Loops over tuples as (variable name, timeseries).

        :param return_as: how to return columns: `xr` as xr.DataArray, `signal` as
            instance of NeuroSignal with the same attributes as the mother signal
        :type return_as: str
        """
        try:
            stacked = self.data.stack({"all": self.dims_not_time})
        except ValueError:
            logging.warning("No dimensions along which to stack...")
            stacked = self.data.expand_dims("all")

        if return_as == "xr":
            yield from stacked.groupby("all")
        elif return_as == "signal":
            for name_coords, column in stacked.groupby("all"):
                if not isinstance(name_coords, (list, tuple)):
                    name_coords = [name_coords]
                name_dict = {k: v for k, v in zip(self.dims_not_time, name_coords)}
                yield name_dict, self.__constructor__(column).__finalize__(self, [f"select {column.name}"])
        else:
            raise ValueError(f"Data type not understood: {return_as}")

    def sel(self, sel_args, inplace=True):
        """
        Subselect part of signal using xarray's `sel`, i.e. selecting by actual
        physical index, hence time in seconds.

        :param sel_args: arguments you'd give to xr.sel(), i.e. slice of times
            you want to select, in seconds as a len=2 list or tuple
        :type sel_args: tuple|list
        :param inplace: whether to do the operation in place or return
        :type inplace: bool
        """
        assert len(sel_args) == 2, "Must provide 2 arguments"
        selected = self.data.sel(time=slice(sel_args[0], sel_args[1]))
        add_steps = [f"select {sel_args[0] or 'x'}:{sel_args[1] or 'x'}s"]
        if inplace:
            self.data = selected
            self.process_steps += add_steps
        else:
            return self.__constructor__(selected).__finalize__(self, add_steps)

    def isel(self, isel_args, inplace=True):
        """
        Subselect part of signal using xarray's `isel`, i.e. selecting by index,
        hence integers.

        :param loc_args: arguments you'd give to xr.isel(), i.e. slice of
            indices you want to select, in seconds as a len=2 list or tuple
        :type loc_args: tuple|list
        :param inplace: whether to do the operation in place or return
        :type inplace: bool
        """
        assert len(isel_args) == 2, "Must provide 2 arguments"
        selected = self.data.isel(time=slice(isel_args[0], isel_args[1]))
        start = isel_args[0] * self.dt if isel_args[0] is not None else "x"
        end = isel_args[1] * self.dt if isel_args[1] is not None else "x"
        add_steps = [f"select {start}:{end}s"]
        if inplace:
            self.data = selected
            self.process_steps += add_steps
        else:
            return self.__constructor__(selected).__finalize__(self, add_steps)

    def rolling(self, roll_over, function=np.mean, dropnans=True, inplace=True):
        """
        Return rolling reduction over signal's time dimension. The window is
        centered around the midpoint.

        :param roll_over: window to use, in seconds
        :type roll_over: float
        :param function: function to use for reduction
        :type function: callable
        :param dropnans: whether to drop NaNs - will shorten time dimension, or
            not
        :type dropnans: bool
        :param inplace: whether to do the operation in place or return
        :type inplace: bool
        """
        assert callable(function)
        rolling = self.data.rolling(time=int(roll_over * self.sampling_frequency), center=True).reduce(function)
        add_steps = [f"rolling {function.__name__} over {roll_over}s"]
        if dropnans:
            rolling = rolling.dropna("time")
            add_steps[0] += "; drop NaNs"
        if inplace:
            self.data = rolling
            self.process_steps += add_steps
        else:
            return self.__constructor__(rolling).__finalize__(self, add_steps)

    def sliding_window(self, length, step=1, window_function="boxcar", lengths_in_seconds=False):
        """
        Return iterator over sliding windows with windowing function applied.
        Each window has length `length` and each is translated by `step` steps.
        For no windowing function use "boxcar". If the last window would have
        the same length as other, it is omitted, i.e. last window does not have
        to end with the final timeseries point!

        :param length: length of the window, can be index or time in seconds,
            see `lengths_in_seconds`
        :type length: int|float
        :param step: how much to translate window in the temporal sense, can be
            index or time in seconds, see `lengths_in_seconds`
        :type step: int|float
        :param window_function: windowing function to use, this is passed to
            `get_window()`; see `scipy.signal.windows.get_window` documentation
        :type window_function: str|tuple|float
        :param lengths_in_seconds: if True, `length` and `step` are interpreted
            in seconds, if False they are indices
        :type lengths_in_seconds: bool
        :yield: generator with windowed Signals
        """
        if lengths_in_seconds:
            length = int(length / self.dt)
            step = int(step / self.dt)
        assert (
            length < self.data.time.shape[0]
        ), f"Length must be smaller than time span of the timeseries: {self.data.time.shape[0]}"
        assert step <= length, "Step cannot be larger than length, some part of timeseries would be omitted!"
        current_idx = 0
        add_steps = f"{str(window_function)} window: "
        windowing_function = get_window(window_function, Nx=length)
        while current_idx <= (self.data.time.shape[0] - length):
            yield self.__constructor__(
                self.data.isel(time=slice(current_idx, current_idx + length)) * windowing_function
            ).__finalize__(self, [add_steps + f"{current_idx}:{current_idx + length}"])
            current_idx += step

    @property
    def shape(self):
        """
        Return shape of the data. Time axis is the first one.
        """
        return self.data.shape

    @property
    def dims_not_time(self):
        """
        Return list of dimensions that are not time.
        """
        return [dim for dim in self.data.dims if dim != "time"]

    @property
    def coords_not_time(self):
        """
        Return dict with all coordinates except time.
        """
        return {k: v.values for k, v in self.data.coords.items() if k != "time"}

    @property
    def start_time(self):
        """
        Return starting time of the signal.
        """
        return self.data.time.values[0]

    @property
    def end_time(self):
        """
        Return ending time of the signal.
        """
        return self.data.time.values[-1]

    @property
    def time(self):
        """
        Return time vector.
        """
        return self.data.time.values

    @property
    def preprocessing_steps(self):
        """
        Return preprocessing steps done on the data.
        """
        return " -> ".join(self.process_steps)

    def pad(self, how_much, in_seconds=False, padding_type="constant", side="both", inplace=True, **kwargs):
        """
        Pad signal by `how_much` on given side of given type.

        :param how_much: how much we should pad, can be time points, or seconds,
            see `in_seconds`
        :type how_much: float|int
        :param in_seconds: whether `how_much` is in seconds, if False, it is
            number of time points
        :type in_seconds: bool
        :param padding_type: how to pad the signal, see `np.pad` documentation
        :type padding_type: str
        :param side: which side to pad - "before", "after", or "both"
        :type side: str
        :param inplace: whether to do the operation in place or return
        :type inplace: bool
        :kwargs: passed to `np.pad`
        """
        if in_seconds:
            how_much = int(np.around(how_much / self.dt))
        if side == "before":
            pad_width = (how_much, 0)
            pad_times = np.arange(-how_much, 0) * self.dt + self.data.time.values[0]
            new_times = np.concatenate([pad_times, self.data.time.values], axis=0)
        elif side == "after":
            pad_width = (0, how_much)
            pad_times = np.arange(1, how_much + 1) * self.dt + self.data.time.values[-1]
            new_times = np.concatenate([self.data.time.values, pad_times], axis=0)
        elif side == "both":
            pad_width = (how_much, how_much)
            pad_before = np.arange(-how_much, 0) * self.dt + self.data.time.values[0]
            pad_after = np.arange(1, how_much + 1) * self.dt + self.data.time.values[-1]
            new_times = np.concatenate([pad_before, self.data.time.values, pad_after], axis=0)
            side += " sides"
        else:
            raise ValueError(f"Unknown padding side: {side}")
        # add padding for other axes than time - zeroes
        pad_width = [(0, 0)] * len(self.dims_not_time) + [pad_width]
        padded = np.pad(self.data.values, pad_width, mode=padding_type, **kwargs)
        # to dataframe
        padded = xr.DataArray(padded, dims=self.data.dims, coords={**self.coords_not_time, "time": new_times})
        add_steps = [f"{how_much * self.dt}s {padding_type} {side} padding"]
        if inplace:
            self.data = padded
            self.process_steps += add_steps
        else:
            return self.__constructor__(padded).__finalize__(self, add_steps)

    def normalize(self, std=False, inplace=True):
        """
        De-mean the timeseries. Optionally also standardise.

        :param std: normalize by std, i.e. to unit variance
        :type std: bool
        :param inplace: whether to do the operation in place or return
        :type inplace: bool
        """

        def norm_func(x, dim):
            demeaned = x - x.mean(dim=dim)
            if std:
                return demeaned / x.std(dim=dim)
            else:
                return demeaned

        normalized = norm_func(self.data, dim="time")
        add_steps = ["normalize", "standardize"] if std else ["normalize"]
        if inplace:
            self.data = normalized
            self.process_steps += add_steps
        else:
            return self.__constructor__(normalized).__finalize__(self, add_steps)

    def resample(self, to_frequency, inplace=True):
        """
        Resample signal to target frequency.

        :param to_frequency: target frequency of the signal, in Hz
        :type to_frequency: float
        :param inplace: whether to do the operation in place or return
        :type inplace: bool
        """
        to_frequency = float(to_frequency)
        try:
            from mne.filter import resample

            resample_func = partial(
                resample, up=to_frequency, down=self.sampling_frequency, npad="auto", axis=-1, pad="edge"
            )
        except ImportError:
            logging.warning("`mne` module not found, falling back to basic scipy's function")

            def resample_func(x):
                return scipy_resample(
                    x,
                    num=int(round((to_frequency / self.sampling_frequency) * self.data.shape[-1])),
                    axis=-1,
                    window="boxcar",
                )

        resampled = resample_func(self.data.values)
        # construct new times
        new_times = (np.arange(resampled.shape[-1], dtype=float) / to_frequency) + self.data.time.values[0]
        # to dataframe
        resampled = xr.DataArray(resampled, dims=self.data.dims, coords={**self.coords_not_time, "time": new_times})
        add_steps = [f"resample to {to_frequency}Hz"]
        if inplace:
            self.data = resampled
            self.sampling_frequency = to_frequency
            self.dt = np.around(np.diff(resampled.time).mean(), 6)
            self.process_steps += add_steps
        else:
            return self.__constructor__(resampled).__finalize__(self, add_steps)

    def hilbert_transform(self, return_as="complex", inplace=True):
        """
        Perform hilbert transform on the signal resulting in analytic signal.

        :param return_as: what to return
            `complex` will compute only analytical signal
            `amplitude` will compute amplitude, hence abs(H(x))
            `phase_wrapped` will compute phase, hence angle(H(x)), in -pi,pi
            `phase_unwrapped` will compute phase in a continuous sense, hence
                monotonic
        :param inplace: whether to do the operation in place or return
        :type inplace: bool
        """
        analytic = hilbert(self.data, axis=-1)
        if return_as == "amplitude":
            analytic = np.abs(analytic)
            add_steps = ["Hilbert - amplitude"]
        elif return_as == "phase_unwrapped":
            analytic = np.unwrap(np.angle(analytic))
            add_steps = ["Hilbert - unwrapped phase"]
        elif return_as == "phase_wrapped":
            analytic = np.angle(analytic)
            add_steps = ["Hilbert - wrapped phase"]
        elif return_as == "complex":
            add_steps = ["Hilbert - complex"]
        else:
            raise ValueError(f"Do not know how to return: {return_as}")

        analytic = xr.DataArray(analytic, dims=self.data.dims, coords=self.data.coords)
        if inplace:
            self.data = analytic
            self.process_steps += add_steps
        else:
            return self.__constructor__(analytic).__finalize__(self, add_steps)

    def detrend(self, segments=None, inplace=True):
        """
        Linearly detrend signal. If segments are given, detrending will be
        performed in each part.

        :param segments: segments for detrending, if None will detrend whole
            signal, given as indices of the time array
        :type segments: list|None
        :param inplace: whether to do the operation in place or return
        :type inplace: bool
        """
        segments = segments or 0
        detrended = detrend(self.data, type="linear", bp=segments, axis=-1)
        detrended = xr.DataArray(detrended, dims=self.data.dims, coords=self.data.coords)
        segments_text = f" with segments: {segments}" if segments != 0 else ""
        add_steps = [f"detrend{segments_text}"]
        if inplace:
            self.data = detrended
            self.process_steps += add_steps
        else:
            return self.__constructor__(detrended).__finalize__(self, add_steps)

    def filter(self, low_freq, high_freq, l_trans_bandwidth="auto", h_trans_bandwidth="auto", inplace=True, **kwargs):
        """
        Filter data. Can be:
            low-pass (low_freq is None, high_freq is not None),
            high-pass (high_freq is None, low_freq is not None),
            band-pass (l_freq < h_freq),
            band-stop (l_freq > h_freq) filter type

        :param low_freq: frequency below which to filter the data
        :type low_freq: float|None
        :param high_freq: frequency above which to filter the data
        :type high_freq: float|None
        :param l_trans_bandwidth: transition band width for low frequency
        :type l_trans_bandwidth: float|str
        :param h_trans_bandwidth: transition band width for high frequency
        :type h_trans_bandwidth: float|str
        :param inplace: whether to do the operation in place or return
        :type inplace: bool
        :**kwargs: possible keywords to `mne.filter.create_filter`:
            `filter_length`="auto",
            `method`="fir",
            `iir_params`=None
            `phase`="zero",
            `fir_window`="hamming",
            `fir_design`="firwin"
        """
        try:
            from mne.filter import filter_data

        except ImportError:
            logging.warning("`mne` module not found, falling back to basic scipy's function")
            filter_data = scipy_iir_filter_data

        filtered = filter_data(
            self.data.values,  # times has to be the last axis
            sfreq=self.sampling_frequency,
            l_freq=low_freq,
            h_freq=high_freq,
            l_trans_bandwidth=l_trans_bandwidth,
            h_trans_bandwidth=h_trans_bandwidth,
            **kwargs,
        )
        add_steps = [f"filter: low {low_freq or 'x'}Hz - high {high_freq or 'x'}Hz"]
        # to dataframe
        filtered = xr.DataArray(filtered, dims=self.data.dims, coords=self.data.coords)
        if inplace:
            self.data = filtered
            self.process_steps += add_steps
        else:
            return self.__constructor__(filtered).__finalize__(self, add_steps)

    def functional_connectivity(self, fc_function=np.corrcoef):
        """
        Compute and return functional connectivity from the data.

        :param fc_function: function which to use for FC computation, should
            take 2D array as space x time and convert it to space x space with
            desired measure
        """
        if len(self.data["space"]) <= 1:
            logging.error("Cannot compute functional connectivity from one timeseries.")
            return None
        if self.data.ndim == 3:
            assert callable(fc_function)
            fcs = []
            for output in self.data["output"]:
                current_slice = self.data.sel({"output": output})
                assert current_slice.ndim == 2
                fcs.append(fc_function(current_slice.values))

            return xr.DataArray(
                np.array(fcs),
                dims=["output", "space", "space"],
                coords={"output": self.data.coords["output"], "space": self.data.coords["space"]},
            )
        if self.data.ndim == 2:
            return xr.DataArray(
                fc_function(self.data.values),
                dims=["space", "space"],
                coords={"space": self.data.coords["space"]},
            )

    def apply(self, func, inplace=True):
        """
        Apply func for each timeseries.

        :param func: function to be applied for each 1D timeseries
        :type func: callable
        :param inplace: whether to do the operation in place or return
        :type inplace: bool
        """
        assert callable(func)
        try:
            # this will work for element-wise function that does not reduces dimensions
            processed = xr.apply_ufunc(func, self.data, input_core_dims=[["time"]], output_core_dims=[["time"]])
            add_steps = [f"apply `{func.__name__}` function over time dim"]
            if inplace:
                self.data = processed
                self.process_steps += add_steps
            else:
                return self.__constructor__(processed).__finalize__(self, add_steps)
        except ValueError:
            # this works for functions that reduce time dimension
            processed = xr.apply_ufunc(func, self.data, input_core_dims=[["time"]])
            logging.warning(
                f"Shape changed after operation! Old shape: {self.shape}, new "
                f"shape: {processed.shape}; Cannot cast to Signal class, "
                "returing as `xr.DataArray`"
            )
            return processed

__constructor__ property

Return constructor, so that each child class would initiate a new instance of the correct class, i.e. first in the method resolution order.

coords_not_time property

Return dict with all coordinates except time.

dims_not_time property

Return list of dimensions that are not time.

end_time property

Return ending time of the signal.

preprocessing_steps property

Return preprocessing steps done on the data.

shape property

Return shape of the data. Time axis is the first one.

start_time property

Return starting time of the signal.

time property

Return time vector.

__eq__(other)

Comparison operator.

Parameters:

Name Type Description Default
other `Signal`

other Signal to compare with

required

Returns:

Type Description
bool

whether two Signals are the same

Source code in neurolib/utils/signal.py
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
def __eq__(self, other):
    """
    Comparison operator.

    :param other: other `Signal` to compare with
    :type other: `Signal`
    :return: whether two `Signals` are the same
    :rtype: bool
    """
    assert isinstance(other, Signal)
    # assert data are the same
    try:
        xr.testing.assert_allclose(self.data, other.data)
        eq = True
    except AssertionError:
        eq = False
    # check attributes, but if not equal, only warn the user
    for attr in self._copy_attributes:
        if getattr(self, attr) != getattr(other, attr):
            logging.warning(f"`{attr}` not equal between signals.")
    return eq

__finalize__(other, add_steps=None)

Copy attributes from other to self. Used when constructing class instance with different data, but same metadata.

Parameters:

Name Type Description Default
other `Signal`

other instance of Signal

required
add_steps list|None

add steps to preprocessing

None
Source code in neurolib/utils/signal.py
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
def __finalize__(self, other, add_steps=None):
    """
    Copy attributes from other to self. Used when constructing class
    instance with different data, but same metadata.

    :param other: other instance of `Signal`
    :type other: `Signal`
    :param add_steps: add steps to preprocessing
    :type add_steps: list|None
    """
    assert isinstance(other, Signal)
    for attr in self._copy_attributes:
        setattr(self, attr, deepcopy(getattr(other, attr)))
    if add_steps is not None:
        self.process_steps += add_steps
    return self

__getitem__(pos)

Get item selects in output dimension.

Source code in neurolib/utils/signal.py
177
178
179
180
181
182
def __getitem__(self, pos):
    """
    Get item selects in output dimension.
    """
    add_steps = [f"select `{pos}` output"]
    return self.__constructor__(self.data.sel(output=pos)).__finalize__(self, add_steps)

__init__(data, time_in_ms=False)

Parameters:

Name Type Description Default
data xr.DataArray

data for the signal, assumes time dimension with time in seconds

required
time_in_ms bool

whether time dimension is in ms

False
Source code in neurolib/utils/signal.py
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
def __init__(self, data, time_in_ms=False):
    """
    :param data: data for the signal, assumes time dimension with time in seconds
    :type data: xr.DataArray
    :param time_in_ms: whether time dimension is in ms
    :type time_in_ms: bool
    """
    assert isinstance(data, xr.DataArray)
    data = deepcopy(data)
    assert "time" in data.dims, "DataArray must have time axis"
    if time_in_ms:
        data["time"] = data["time"] / 1000.0
    data["time"] = np.around(data["time"], 6)
    self.data = data
    # assert time dimension is last
    self.data = self.data.transpose(*(self.dims_not_time + ["time"]))
    # compute dt and sampling frequency
    self.dt = np.around(np.diff(data.time).mean(), 6)
    self.sampling_frequency = 1.0 / self.dt
    self.process_steps = [f"raw {self.signal_type} signal: {self.start_time}--{self.end_time}s"]

__repr__()

Representation.

Source code in neurolib/utils/signal.py
149
150
151
152
153
def __repr__(self):
    """
    Representation.
    """
    return self.__str__()

__str__()

String representation.

Source code in neurolib/utils/signal.py
137
138
139
140
141
142
143
144
145
146
147
def __str__(self):
    """
    String representation.
    """
    return (
        f"{self.name} representing {self.signal_type} signal with unit of "
        f"{self.unit} with user-provided description: `{self.description}`"
        f". Shape of the signal is {self.shape} with dimensions "
        f"{self.data.dims}. Signal starts at {self.start_time} and ends at "
        f"{self.end_time}."
    )

apply(func, inplace=True)

Apply func for each timeseries.

Parameters:

Name Type Description Default
func callable

function to be applied for each 1D timeseries

required
inplace bool

whether to do the operation in place or return

True
Source code in neurolib/utils/signal.py
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
def apply(self, func, inplace=True):
    """
    Apply func for each timeseries.

    :param func: function to be applied for each 1D timeseries
    :type func: callable
    :param inplace: whether to do the operation in place or return
    :type inplace: bool
    """
    assert callable(func)
    try:
        # this will work for element-wise function that does not reduces dimensions
        processed = xr.apply_ufunc(func, self.data, input_core_dims=[["time"]], output_core_dims=[["time"]])
        add_steps = [f"apply `{func.__name__}` function over time dim"]
        if inplace:
            self.data = processed
            self.process_steps += add_steps
        else:
            return self.__constructor__(processed).__finalize__(self, add_steps)
    except ValueError:
        # this works for functions that reduce time dimension
        processed = xr.apply_ufunc(func, self.data, input_core_dims=[["time"]])
        logging.warning(
            f"Shape changed after operation! Old shape: {self.shape}, new "
            f"shape: {processed.shape}; Cannot cast to Signal class, "
            "returing as `xr.DataArray`"
        )
        return processed

detrend(segments=None, inplace=True)

Linearly detrend signal. If segments are given, detrending will be performed in each part.

Parameters:

Name Type Description Default
segments list|None

segments for detrending, if None will detrend whole signal, given as indices of the time array

None
inplace bool

whether to do the operation in place or return

True
Source code in neurolib/utils/signal.py
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
def detrend(self, segments=None, inplace=True):
    """
    Linearly detrend signal. If segments are given, detrending will be
    performed in each part.

    :param segments: segments for detrending, if None will detrend whole
        signal, given as indices of the time array
    :type segments: list|None
    :param inplace: whether to do the operation in place or return
    :type inplace: bool
    """
    segments = segments or 0
    detrended = detrend(self.data, type="linear", bp=segments, axis=-1)
    detrended = xr.DataArray(detrended, dims=self.data.dims, coords=self.data.coords)
    segments_text = f" with segments: {segments}" if segments != 0 else ""
    add_steps = [f"detrend{segments_text}"]
    if inplace:
        self.data = detrended
        self.process_steps += add_steps
    else:
        return self.__constructor__(detrended).__finalize__(self, add_steps)

filter(low_freq, high_freq, l_trans_bandwidth='auto', h_trans_bandwidth='auto', inplace=True, **kwargs)

Filter data. Can be: low-pass (low_freq is None, high_freq is not None), high-pass (high_freq is None, low_freq is not None), band-pass (l_freq < h_freq), band-stop (l_freq > h_freq) filter type

:**kwargs: possible keywords to mne.filter.create_filter: filter_length="auto", method="fir", iir_params=None phase="zero", fir_window="hamming", fir_design="firwin"

Parameters:

Name Type Description Default
low_freq float|None

frequency below which to filter the data

required
high_freq float|None

frequency above which to filter the data

required
l_trans_bandwidth float|str

transition band width for low frequency

'auto'
h_trans_bandwidth float|str

transition band width for high frequency

'auto'
inplace bool

whether to do the operation in place or return

True
Source code in neurolib/utils/signal.py
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
636
637
638
def filter(self, low_freq, high_freq, l_trans_bandwidth="auto", h_trans_bandwidth="auto", inplace=True, **kwargs):
    """
    Filter data. Can be:
        low-pass (low_freq is None, high_freq is not None),
        high-pass (high_freq is None, low_freq is not None),
        band-pass (l_freq < h_freq),
        band-stop (l_freq > h_freq) filter type

    :param low_freq: frequency below which to filter the data
    :type low_freq: float|None
    :param high_freq: frequency above which to filter the data
    :type high_freq: float|None
    :param l_trans_bandwidth: transition band width for low frequency
    :type l_trans_bandwidth: float|str
    :param h_trans_bandwidth: transition band width for high frequency
    :type h_trans_bandwidth: float|str
    :param inplace: whether to do the operation in place or return
    :type inplace: bool
    :**kwargs: possible keywords to `mne.filter.create_filter`:
        `filter_length`="auto",
        `method`="fir",
        `iir_params`=None
        `phase`="zero",
        `fir_window`="hamming",
        `fir_design`="firwin"
    """
    try:
        from mne.filter import filter_data

    except ImportError:
        logging.warning("`mne` module not found, falling back to basic scipy's function")
        filter_data = scipy_iir_filter_data

    filtered = filter_data(
        self.data.values,  # times has to be the last axis
        sfreq=self.sampling_frequency,
        l_freq=low_freq,
        h_freq=high_freq,
        l_trans_bandwidth=l_trans_bandwidth,
        h_trans_bandwidth=h_trans_bandwidth,
        **kwargs,
    )
    add_steps = [f"filter: low {low_freq or 'x'}Hz - high {high_freq or 'x'}Hz"]
    # to dataframe
    filtered = xr.DataArray(filtered, dims=self.data.dims, coords=self.data.coords)
    if inplace:
        self.data = filtered
        self.process_steps += add_steps
    else:
        return self.__constructor__(filtered).__finalize__(self, add_steps)

from_file(filename) classmethod

Load signal from saved file.

Parameters:

Name Type Description Default
filename str

filename for the Signal

required
Source code in neurolib/utils/signal.py
 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
@classmethod
def from_file(cls, filename):
    """
    Load signal from saved file.

    :param filename: filename for the Signal
    :type filename: str
    """
    if not filename.endswith(NC_EXT):
        filename += NC_EXT
    # load NC file
    xarray = xr.load_dataarray(filename)
    # init class
    signal = cls(xarray)
    # if nc file has attributes, copy them to signal class
    if xarray.attrs:
        process_steps = []
        for k, v in xarray.attrs.items():
            if cls.PROCESS_STEPS_KEY in k:
                idx = int(k[len(cls.PROCESS_STEPS_KEY) + 1 :])
                process_steps.insert(idx, v)
            else:
                setattr(signal, k, v)
    else:
        logging.warning("No metadata found, setting empty...")
        process_steps = [f"raw {signal.signal_type} signal: {signal.start_time}--" f"{signal.end_time}s"]
    setattr(signal, cls.PROCESS_STEPS_KEY, process_steps)
    return signal

from_model_output(model, group='', time_in_ms=True) classmethod

Initial Signal from modelling output.

Source code in neurolib/utils/signal.py
79
80
81
82
83
84
85
@classmethod
def from_model_output(cls, model, group="", time_in_ms=True):
    """
    Initial Signal from modelling output.
    """
    assert isinstance(model, Model)
    return cls(model.xr(group=group), time_in_ms=time_in_ms)

functional_connectivity(fc_function=np.corrcoef)

Compute and return functional connectivity from the data.

Parameters:

Name Type Description Default
fc_function

function which to use for FC computation, should take 2D array as space x time and convert it to space x space with desired measure

np.corrcoef
Source code in neurolib/utils/signal.py
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
def functional_connectivity(self, fc_function=np.corrcoef):
    """
    Compute and return functional connectivity from the data.

    :param fc_function: function which to use for FC computation, should
        take 2D array as space x time and convert it to space x space with
        desired measure
    """
    if len(self.data["space"]) <= 1:
        logging.error("Cannot compute functional connectivity from one timeseries.")
        return None
    if self.data.ndim == 3:
        assert callable(fc_function)
        fcs = []
        for output in self.data["output"]:
            current_slice = self.data.sel({"output": output})
            assert current_slice.ndim == 2
            fcs.append(fc_function(current_slice.values))

        return xr.DataArray(
            np.array(fcs),
            dims=["output", "space", "space"],
            coords={"output": self.data.coords["output"], "space": self.data.coords["space"]},
        )
    if self.data.ndim == 2:
        return xr.DataArray(
            fc_function(self.data.values),
            dims=["space", "space"],
            coords={"space": self.data.coords["space"]},
        )

hilbert_transform(return_as='complex', inplace=True)

Perform hilbert transform on the signal resulting in analytic signal.

Parameters:

Name Type Description Default
return_as

what to return complex will compute only analytical signal amplitude will compute amplitude, hence abs(H(x)) phase_wrapped will compute phase, hence angle(H(x)), in -pi,pi phase_unwrapped will compute phase in a continuous sense, hence monotonic

'complex'
inplace bool

whether to do the operation in place or return

True
Source code in neurolib/utils/signal.py
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
def hilbert_transform(self, return_as="complex", inplace=True):
    """
    Perform hilbert transform on the signal resulting in analytic signal.

    :param return_as: what to return
        `complex` will compute only analytical signal
        `amplitude` will compute amplitude, hence abs(H(x))
        `phase_wrapped` will compute phase, hence angle(H(x)), in -pi,pi
        `phase_unwrapped` will compute phase in a continuous sense, hence
            monotonic
    :param inplace: whether to do the operation in place or return
    :type inplace: bool
    """
    analytic = hilbert(self.data, axis=-1)
    if return_as == "amplitude":
        analytic = np.abs(analytic)
        add_steps = ["Hilbert - amplitude"]
    elif return_as == "phase_unwrapped":
        analytic = np.unwrap(np.angle(analytic))
        add_steps = ["Hilbert - unwrapped phase"]
    elif return_as == "phase_wrapped":
        analytic = np.angle(analytic)
        add_steps = ["Hilbert - wrapped phase"]
    elif return_as == "complex":
        add_steps = ["Hilbert - complex"]
    else:
        raise ValueError(f"Do not know how to return: {return_as}")

    analytic = xr.DataArray(analytic, dims=self.data.dims, coords=self.data.coords)
    if inplace:
        self.data = analytic
        self.process_steps += add_steps
    else:
        return self.__constructor__(analytic).__finalize__(self, add_steps)

isel(isel_args, inplace=True)

Subselect part of signal using xarray's isel, i.e. selecting by index, hence integers.

Parameters:

Name Type Description Default
loc_args tuple|list

arguments you'd give to xr.isel(), i.e. slice of indices you want to select, in seconds as a len=2 list or tuple

required
inplace bool

whether to do the operation in place or return

True
Source code in neurolib/utils/signal.py
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
def isel(self, isel_args, inplace=True):
    """
    Subselect part of signal using xarray's `isel`, i.e. selecting by index,
    hence integers.

    :param loc_args: arguments you'd give to xr.isel(), i.e. slice of
        indices you want to select, in seconds as a len=2 list or tuple
    :type loc_args: tuple|list
    :param inplace: whether to do the operation in place or return
    :type inplace: bool
    """
    assert len(isel_args) == 2, "Must provide 2 arguments"
    selected = self.data.isel(time=slice(isel_args[0], isel_args[1]))
    start = isel_args[0] * self.dt if isel_args[0] is not None else "x"
    end = isel_args[1] * self.dt if isel_args[1] is not None else "x"
    add_steps = [f"select {start}:{end}s"]
    if inplace:
        self.data = selected
        self.process_steps += add_steps
    else:
        return self.__constructor__(selected).__finalize__(self, add_steps)

iterate(return_as='signal')

Return iterator over columns, so univariate measures can be computed per column. Loops over tuples as (variable name, timeseries).

Parameters:

Name Type Description Default
return_as str

how to return columns: xr as xr.DataArray, signal as instance of NeuroSignal with the same attributes as the mother signal

'signal'
Source code in neurolib/utils/signal.py
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
def iterate(self, return_as="signal"):
    """
    Return iterator over columns, so univariate measures can be computed
    per column. Loops over tuples as (variable name, timeseries).

    :param return_as: how to return columns: `xr` as xr.DataArray, `signal` as
        instance of NeuroSignal with the same attributes as the mother signal
    :type return_as: str
    """
    try:
        stacked = self.data.stack({"all": self.dims_not_time})
    except ValueError:
        logging.warning("No dimensions along which to stack...")
        stacked = self.data.expand_dims("all")

    if return_as == "xr":
        yield from stacked.groupby("all")
    elif return_as == "signal":
        for name_coords, column in stacked.groupby("all"):
            if not isinstance(name_coords, (list, tuple)):
                name_coords = [name_coords]
            name_dict = {k: v for k, v in zip(self.dims_not_time, name_coords)}
            yield name_dict, self.__constructor__(column).__finalize__(self, [f"select {column.name}"])
    else:
        raise ValueError(f"Data type not understood: {return_as}")

normalize(std=False, inplace=True)

De-mean the timeseries. Optionally also standardise.

Parameters:

Name Type Description Default
std bool

normalize by std, i.e. to unit variance

False
inplace bool

whether to do the operation in place or return

True
Source code in neurolib/utils/signal.py
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
def normalize(self, std=False, inplace=True):
    """
    De-mean the timeseries. Optionally also standardise.

    :param std: normalize by std, i.e. to unit variance
    :type std: bool
    :param inplace: whether to do the operation in place or return
    :type inplace: bool
    """

    def norm_func(x, dim):
        demeaned = x - x.mean(dim=dim)
        if std:
            return demeaned / x.std(dim=dim)
        else:
            return demeaned

    normalized = norm_func(self.data, dim="time")
    add_steps = ["normalize", "standardize"] if std else ["normalize"]
    if inplace:
        self.data = normalized
        self.process_steps += add_steps
    else:
        return self.__constructor__(normalized).__finalize__(self, add_steps)

pad(how_much, in_seconds=False, padding_type='constant', side='both', inplace=True, **kwargs)

Pad signal by how_much on given side of given type.

:kwargs: passed to np.pad

Parameters:

Name Type Description Default
how_much float|int

how much we should pad, can be time points, or seconds, see in_seconds

required
in_seconds bool

whether how_much is in seconds, if False, it is number of time points

False
padding_type str

how to pad the signal, see np.pad documentation

'constant'
side str

which side to pad - "before", "after", or "both"

'both'
inplace bool

whether to do the operation in place or return

True
Source code in neurolib/utils/signal.py
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
def pad(self, how_much, in_seconds=False, padding_type="constant", side="both", inplace=True, **kwargs):
    """
    Pad signal by `how_much` on given side of given type.

    :param how_much: how much we should pad, can be time points, or seconds,
        see `in_seconds`
    :type how_much: float|int
    :param in_seconds: whether `how_much` is in seconds, if False, it is
        number of time points
    :type in_seconds: bool
    :param padding_type: how to pad the signal, see `np.pad` documentation
    :type padding_type: str
    :param side: which side to pad - "before", "after", or "both"
    :type side: str
    :param inplace: whether to do the operation in place or return
    :type inplace: bool
    :kwargs: passed to `np.pad`
    """
    if in_seconds:
        how_much = int(np.around(how_much / self.dt))
    if side == "before":
        pad_width = (how_much, 0)
        pad_times = np.arange(-how_much, 0) * self.dt + self.data.time.values[0]
        new_times = np.concatenate([pad_times, self.data.time.values], axis=0)
    elif side == "after":
        pad_width = (0, how_much)
        pad_times = np.arange(1, how_much + 1) * self.dt + self.data.time.values[-1]
        new_times = np.concatenate([self.data.time.values, pad_times], axis=0)
    elif side == "both":
        pad_width = (how_much, how_much)
        pad_before = np.arange(-how_much, 0) * self.dt + self.data.time.values[0]
        pad_after = np.arange(1, how_much + 1) * self.dt + self.data.time.values[-1]
        new_times = np.concatenate([pad_before, self.data.time.values, pad_after], axis=0)
        side += " sides"
    else:
        raise ValueError(f"Unknown padding side: {side}")
    # add padding for other axes than time - zeroes
    pad_width = [(0, 0)] * len(self.dims_not_time) + [pad_width]
    padded = np.pad(self.data.values, pad_width, mode=padding_type, **kwargs)
    # to dataframe
    padded = xr.DataArray(padded, dims=self.data.dims, coords={**self.coords_not_time, "time": new_times})
    add_steps = [f"{how_much * self.dt}s {padding_type} {side} padding"]
    if inplace:
        self.data = padded
        self.process_steps += add_steps
    else:
        return self.__constructor__(padded).__finalize__(self, add_steps)

resample(to_frequency, inplace=True)

Resample signal to target frequency.

Parameters:

Name Type Description Default
to_frequency float

target frequency of the signal, in Hz

required
inplace bool

whether to do the operation in place or return

True
Source code in neurolib/utils/signal.py
491
492
493
494
495
496
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
def resample(self, to_frequency, inplace=True):
    """
    Resample signal to target frequency.

    :param to_frequency: target frequency of the signal, in Hz
    :type to_frequency: float
    :param inplace: whether to do the operation in place or return
    :type inplace: bool
    """
    to_frequency = float(to_frequency)
    try:
        from mne.filter import resample

        resample_func = partial(
            resample, up=to_frequency, down=self.sampling_frequency, npad="auto", axis=-1, pad="edge"
        )
    except ImportError:
        logging.warning("`mne` module not found, falling back to basic scipy's function")

        def resample_func(x):
            return scipy_resample(
                x,
                num=int(round((to_frequency / self.sampling_frequency) * self.data.shape[-1])),
                axis=-1,
                window="boxcar",
            )

    resampled = resample_func(self.data.values)
    # construct new times
    new_times = (np.arange(resampled.shape[-1], dtype=float) / to_frequency) + self.data.time.values[0]
    # to dataframe
    resampled = xr.DataArray(resampled, dims=self.data.dims, coords={**self.coords_not_time, "time": new_times})
    add_steps = [f"resample to {to_frequency}Hz"]
    if inplace:
        self.data = resampled
        self.sampling_frequency = to_frequency
        self.dt = np.around(np.diff(resampled.time).mean(), 6)
        self.process_steps += add_steps
    else:
        return self.__constructor__(resampled).__finalize__(self, add_steps)

rolling(roll_over, function=np.mean, dropnans=True, inplace=True)

Return rolling reduction over signal's time dimension. The window is centered around the midpoint.

Parameters:

Name Type Description Default
roll_over float

window to use, in seconds

required
function callable

function to use for reduction

np.mean
dropnans bool

whether to drop NaNs - will shorten time dimension, or not

True
inplace bool

whether to do the operation in place or return

True
Source code in neurolib/utils/signal.py
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
def rolling(self, roll_over, function=np.mean, dropnans=True, inplace=True):
    """
    Return rolling reduction over signal's time dimension. The window is
    centered around the midpoint.

    :param roll_over: window to use, in seconds
    :type roll_over: float
    :param function: function to use for reduction
    :type function: callable
    :param dropnans: whether to drop NaNs - will shorten time dimension, or
        not
    :type dropnans: bool
    :param inplace: whether to do the operation in place or return
    :type inplace: bool
    """
    assert callable(function)
    rolling = self.data.rolling(time=int(roll_over * self.sampling_frequency), center=True).reduce(function)
    add_steps = [f"rolling {function.__name__} over {roll_over}s"]
    if dropnans:
        rolling = rolling.dropna("time")
        add_steps[0] += "; drop NaNs"
    if inplace:
        self.data = rolling
        self.process_steps += add_steps
    else:
        return self.__constructor__(rolling).__finalize__(self, add_steps)

save(filename)

Save signal.

Parameters:

Name Type Description Default
filename str

filename to save, currently saves to netCDF file, which is natively supported by xarray

required
Source code in neurolib/utils/signal.py
224
225
226
227
228
229
230
231
232
233
234
def save(self, filename):
    """
    Save signal.

    :param filename: filename to save, currently saves to netCDF file, which is natively supported by xarray
    :type filename: str
    """
    self._write_attrs_to_xr()
    if not filename.endswith(NC_EXT):
        filename += NC_EXT
    self.data.to_netcdf(filename)

sel(sel_args, inplace=True)

Subselect part of signal using xarray's sel, i.e. selecting by actual physical index, hence time in seconds.

Parameters:

Name Type Description Default
sel_args tuple|list

arguments you'd give to xr.sel(), i.e. slice of times you want to select, in seconds as a len=2 list or tuple

required
inplace bool

whether to do the operation in place or return

True
Source code in neurolib/utils/signal.py
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
def sel(self, sel_args, inplace=True):
    """
    Subselect part of signal using xarray's `sel`, i.e. selecting by actual
    physical index, hence time in seconds.

    :param sel_args: arguments you'd give to xr.sel(), i.e. slice of times
        you want to select, in seconds as a len=2 list or tuple
    :type sel_args: tuple|list
    :param inplace: whether to do the operation in place or return
    :type inplace: bool
    """
    assert len(sel_args) == 2, "Must provide 2 arguments"
    selected = self.data.sel(time=slice(sel_args[0], sel_args[1]))
    add_steps = [f"select {sel_args[0] or 'x'}:{sel_args[1] or 'x'}s"]
    if inplace:
        self.data = selected
        self.process_steps += add_steps
    else:
        return self.__constructor__(selected).__finalize__(self, add_steps)

sliding_window(length, step=1, window_function='boxcar', lengths_in_seconds=False)

Return iterator over sliding windows with windowing function applied. Each window has length length and each is translated by step steps. For no windowing function use "boxcar". If the last window would have the same length as other, it is omitted, i.e. last window does not have to end with the final timeseries point!

:yield: generator with windowed Signals

Parameters:

Name Type Description Default
length int|float

length of the window, can be index or time in seconds, see lengths_in_seconds

required
step int|float

how much to translate window in the temporal sense, can be index or time in seconds, see lengths_in_seconds

1
window_function str|tuple|float

windowing function to use, this is passed to get_window(); see scipy.signal.windows.get_window documentation

'boxcar'
lengths_in_seconds bool

if True, length and step are interpreted in seconds, if False they are indices

False
Source code in neurolib/utils/signal.py
331
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
def sliding_window(self, length, step=1, window_function="boxcar", lengths_in_seconds=False):
    """
    Return iterator over sliding windows with windowing function applied.
    Each window has length `length` and each is translated by `step` steps.
    For no windowing function use "boxcar". If the last window would have
    the same length as other, it is omitted, i.e. last window does not have
    to end with the final timeseries point!

    :param length: length of the window, can be index or time in seconds,
        see `lengths_in_seconds`
    :type length: int|float
    :param step: how much to translate window in the temporal sense, can be
        index or time in seconds, see `lengths_in_seconds`
    :type step: int|float
    :param window_function: windowing function to use, this is passed to
        `get_window()`; see `scipy.signal.windows.get_window` documentation
    :type window_function: str|tuple|float
    :param lengths_in_seconds: if True, `length` and `step` are interpreted
        in seconds, if False they are indices
    :type lengths_in_seconds: bool
    :yield: generator with windowed Signals
    """
    if lengths_in_seconds:
        length = int(length / self.dt)
        step = int(step / self.dt)
    assert (
        length < self.data.time.shape[0]
    ), f"Length must be smaller than time span of the timeseries: {self.data.time.shape[0]}"
    assert step <= length, "Step cannot be larger than length, some part of timeseries would be omitted!"
    current_idx = 0
    add_steps = f"{str(window_function)} window: "
    windowing_function = get_window(window_function, Nx=length)
    while current_idx <= (self.data.time.shape[0] - length):
        yield self.__constructor__(
            self.data.isel(time=slice(current_idx, current_idx + length)) * windowing_function
        ).__finalize__(self, [add_steps + f"{current_idx}:{current_idx + length}"])
        current_idx += step