toad.shifts¶
Shifts module for TOAD.
This module provides functionality for detecting abrupt shifts in climate data. The main function compute_shifts takes a dataset and a method for detecting abrupt shifts and returns the detected shifts as an xarray object.
Currently implemented methods: - ASDETECT: Implementation of the [Boulton+Lenton2019]
Functions
|
Apply an abrupt shift detection algorithm to a dataset along the specified temporal dimension. |
- class toad.shifts.ASDETECT(lmin=5, lmax=None, timescale=None, segmentation='two_sided', ignore_nan_warnings=False)¶
Bases:
ShiftsMethodDetect abrupt shifts in a time series using gradient-based analysis related to [Boulton+Lenton2019].
- The algorithm works by:
Dividing the time series into segments of size l.
Performing linear regression within each segment to calculate gradients.
Identifying significant gradients exceeding ±3 Median Absolute Deviations (MAD) from the median gradient.
Updating a detection array by adding +1 for significant positive gradients and -1 for significant negative gradients in each segment.
Iterating over multiple window sizes (l), updating the detection array at each step.
Normalizing the detection array according to the selected segmentation mode.
Two segmentation modes are available:
- “two_sided” mode (default, recommended):
Uses forward and backward segmentation to reduce bias and improve coverage.
Segments overlap, providing better temporal coverage of the time series.
Normalizes by dividing by a counter array that tracks how many segments cover each time point.
Applies edge correction by downweighting signals at boundaries where fewer segments overlap, reducing unreliable detections at the start and end of the time series.
Generally produces smoother and more reliable results, especially near the edges.
- “original” mode (following [Boulton+Lenton2019]):
Uses centered, non-overlapping segments that are evenly distributed within the time series.
Normalizes by dividing by the total number of window sizes used (lmax - lmin + 1).
Simpler approach that matches the original algorithm implementation.
Note: ASDETECT does not work with NaN values so it will return a detection time series of all zeros if the input time series contains NaN values.
- Parameters:
lmin (int) – The minimum segment length (in time steps) used for detection. Controls the shortest resolvable timescale and must be at least 3. Defaults to 5.
lmax (int | None) – Optional maximum segment length (in time steps) used for detection. To resolve multiple shifts within the same grid cell, this should be chosen smaller than the minimum temporal separation between shifts of interest. If not specified, it defaults to one third of the length of the time dimension.
timescale (tuple[float | None, float | None] | None) –
Optional tuple specifying the minimum and maximum time window sizes for shift detection in the same physical units as the input time axis (e.g. years, days). When provided, this is converted internally to integer segment lengths (lmin, lmax) and overrides the explicit lmin/lmax settings. Either bound may be set to None to fall back to the defaults (5 time steps for the minimum, one third of the series length for the maximum).
Important: When the input dataset uses cftime format, the timescale parameter must be specified in seconds, because cftime values are internally converted to seconds since the first time point. For example, to detect shifts at timescales of 1-2 years with cftime data, use
timescale=(365.25*24*3600, 2*365.25*24*3600)(seconds). For numeric time data (not cftime), use the same units as your time axis.Note: The timescale conversion requires regular temporal spacing. If your time series has irregular spacing, a ValueError will be raised. In this case, use lmin/lmax parameters directly instead of timescale, or resample your data to regular spacing.
segmentation (Literal['two_sided', 'original'] | str) – Segmentation method to use. “two_sided” (recommended) applies forward and backward segmentation with counter-based normalisation and edge correction, reducing positional bias at roughly twice the computational cost of the original algorithm. “original” reproduces the centred segmentation of the original ASDETECT implementation for backward compatibility. Defaults to “two_sided”.
ignore_nan_warnings (bool) – (Optional) If True, timeseries containing NaN values will be ignored, i.e. a detection time series of all zeros will be returned. If False, an error will be raised.
- LMIN_MIN = 5¶
- fit_predict(values_1d, times_1d)¶
Compute the detection time series for each grid cell in the 3D data array.
- Parameters:
values_1d (ndarray[tuple[int, ...], dtype[float64]]) – 1D array of values
times_1d (ndarray[tuple[int, ...], dtype[float64]]) – 1D array of times
- Returns:
1 indicates that all tested segment lengths detected a significant positive gradient (i.e. exceeding 3 MAD of the median gradient),
-1 indicates that all tested segment lengths detected a significant negative gradient.
Values between -1 and 1 indicate the proportion of segment lengths detecting a significant gradient at that time point.
- Return type:
A 1D array of the same length as values_1d, where each value represents the abrupt shift score for a grid cell at a specific time. The score ranges from -1 to 1
- pre_validation(data_array, td)¶
Validate regular temporal spacing and convert timescale to lmin/lmax if needed.
This method is called once in compute_shifts() before applying the method to all grid cells. It always validates that the time series has regular spacing (required for the algorithm to work correctly) and converts timescale parameters to lmin/lmax if timescale is provided.
- Parameters:
data_array (xr.DataArray) – The masked data array that will be processed
td (TOAD) – The TOAD object containing the dataset and metadata
- Raises:
ValueError – If time spacing is irregular
- Return type:
None
- class toad.shifts.ShiftsMethod¶
Bases:
ABC- abstractmethod fit_predict(values_1d, times_1d)¶
Apply the shifts detection method.
- Parameters:
values_1d (ndarray) – Input values.
times_1d (ndarray) – Input times.
- Returns:
- A detection time series with the same length as the input,
where each value indicates the presence or magnitude of a detected shift.
- Return type:
np.ndarray
- pre_validation(data_array, td)¶
Optional validation method that runs once before applying the method to all grid cells.
This method is called once in compute_shifts() before xr.apply_ufunc() processes all grid cells. Override this method in subclasses to implement method-specific validations (e.g., checking for regular temporal spacing, validating parameters, converting timescale parameters, etc.).
- Parameters:
data_array (xr.DataArray) – The masked data array that will be processed
td (TOAD) – The TOAD object containing the dataset and metadata
- Raises:
ValueError – If validation fails
- Return type:
None
- toad.shifts.compute_shifts(td, var, method, output_label_suffix='', overwrite=False)¶
Apply an abrupt shift detection algorithm to a dataset along the specified temporal dimension.
- Parameters:
td (TOAD) – TOAD object
var (str) – Name of the variable in the dataset to analyze for abrupt shifts
method (ShiftsMethod) – The abrupt shift detection algorithm to use. Choose from predefined method objects in toad.shifts or create your own following the base class in toad.shifts.methods.base
output_label_suffix (str) – A suffix to add to the output label. Defaults to “”.
overwrite (bool) – Whether to overwrite existing variable. Defaults to False.
- Returns:
- If merge_input is True, returns an xarray.Dataset containing
the original data and the detected shifts. If merge_input is False, returns an xarray.DataArray containing the detected shifts.
- Return type:
Union[xr.Dataset, xr.DataArray]
- Raises:
ValueError – If data is invalid or required parameters are missing
Modules
Shift detection methods are exposed in shifts/__init__.py |