Bias adjustment and downscaling algorithms

xarray data structures allow for relatively straightforward implementations of simple bias-adjustment and downscaling algorithms documented in Adjustment methods. Each algorithm is split into train and adjust components. The train function will compare two DataArrays x and y, and create a dataset storing the transfer information allowing to go from x to y. This dataset, stored in the adjustment object, can then be used by the adjust method to apply this information to x. x could be the same DataArray used for training, or another DataArray with similar characteristics.

For example, given a daily time series of observations ref, a model simulation over the observational period hist and a model simulation over a future period sim, we would apply a bias-adjustment method such as detrended quantile mapping (DQM) as:

from xclim import sdba
dqm = sdba.adjustment.DetrendedQuantileMapping.train(ref, hist)
scen = dqm.adjust(sim)

Most method can either be applied additively or multiplicatively. Also, most methods can be applied independently on different time groupings (monthly, seasonally) or according to the day of the year and a rolling window width.

When transfer factors are applied in adjustment, they can be interpolated according to the time grouping. This helps avoid discontinuities in adjustment factors at the beginning of each season or month and is computationaly cheaper than computing adjustment factors for each day of the year. (Currently only implemented for monthly grouping)

Application in multivariate settings

When applying univariate adjustment methods to multiple variables, some strategies are recommended to avoid introducing unrealistic artifacts in adjusted outputs.

Minimum and maximum temperature

When adjusting both minimum and maximum temperature, adjustment factors sometimes yield minimum temperatures larger than the maximum temperature on the same day, which of course, is non-sensical. One way to avoid this is to first adjust maximum temperature using an additive adjustment, then adjust the diurnal temperature range (DTR) using a multiplicative adjustment, and then determine minimum temperature by subtracting DTR from the maximum temperature ([Thrasher], [Agbazo])

Relative and specific humidity

When adjusting both relative and specific humidity, we want to preserve the relationship between both. To do this, [Grenier] suggests to first adjust the relative humidity using a multiplicative factor, ensure values are within 0-100%, then apply an additive adjustment factor to the surface pressure before estimating the specific humidity from thermodynamic relationships.

Radiation and precipitation

In theory, short wave radiation should be capped when precipitation is not zero, but there is as of yet no mechanism proposed to do that, see [Hoffman].

References

Agbazo

Agbazo, M. N., & Grenier, P. (2019). Characterizing and avoiding physical inconsistency generated by the application of univariate quantile mapping on daily minimum and maximum temperatures over Hudson Bay. International Journal of Climatology, joc.6432. https://doi.org/10.1002/joc.6432

Grenier

Grenier, P. (2018). Two Types of Physical Inconsistency to Avoid with Univariate Quantile Mapping: A Case Study over North America Concerning Relative Humidity and Its Parent Variables. Journal of Applied Meteorology and Climatology, 57(2), 347–364. https://doi.org/10.1175/JAMC-D-17-0177.1

Hoffman

Hoffmann, H., & Rath, T. (2012). Meteorologically consistent bias correction of climate time series for agricultural models. Theoretical and Applied Climatology, 110(1–2), 129–141. https://doi.org/10.1007/s00704-012-0618-x

Thrasher

Thrasher, B., Maurer, E. P., McKellar, C., & Duffy, P. B. (2012). Technical Note: Bias correcting climate model simulated daily temperature extremes with quantile mapping. Hydrology and Earth System Sciences, 16(9), 3309–3314. https://doi.org/10.5194/hess-16-3309-2012

SDBA’s user API

Adjustment methods

class xclim.sdba.adjustment.DetrendedQuantileMapping(*args, _trained=False, **kwargs)[source]

Bases: xclim.sdba.adjustment.TrainAdjust

Detrended Quantile Mapping bias-adjustment.

The algorithm follows these steps, 1-3 being the ‘train’ and 4-6, the ‘adjust’ steps.

  1. A scaling factor that would make the mean of hist match the mean of ref is computed.

  2. ref and hist are normalized by removing the “dayofyear” mean.

  3. Adjustment factors are computed between the quantiles of the normalized ref and hist.

  4. sim is corrected by the scaling factor, and either normalized by “dayofyear” and detrended group-wise or directly detrended per “dayofyear”, using a linear fit (modifiable).

  5. Values of detrended sim are matched to the corresponding quantiles of normalized hist and corrected accordingly.

  6. The trend is put back on the result.

\[F^{-1}_{ref}\left\{F_{hist}\left[\frac{\overline{hist}\cdot sim}{\overline{sim}}\right]\right\}\frac{\overline{sim}}{\overline{hist}}\]

where \(F\) is the cumulative distribution function (CDF) and \(\overline{xyz}\) is the linear trend of the data. This equation is valid for multiplicative adjustment. Based on the DQM method of [Cannon2015].

Parameters
  • Train step

  • nquantiles (int or 1d array of floats) – The number of quantiles to use. Two endpoints at 1e-6 and 1 - 1e-6 will be added. An array of quantiles [0, 1] can also be passed. Defaults to 20 quantiles (+2 endpoints).

  • kind ({‘+’, ‘’}*) – The adjustment kind, either additive or multiplicative. Defaults to “+”.

  • group (Union[str, Grouper]) – The grouping information. See xclim.sdba.base.Grouper for details. Default is “time”, meaning an single adjustment group along dimension “time”.

  • Adjust step

  • interp ({‘nearest’, ‘linear’, ‘cubic’}) – The interpolation method to use when interpolating the adjustment factors. Defaults to “nearest”.

  • detrend (int or BaseDetrend instance) – The method to use when detrending. If an int is passed, it is understood as a PolyDetrend (polynomial detrending) degree. Defaults to 1 (linear detrending)

  • extrapolation ({‘constant’, ‘nan’}) – The type of extrapolation to use. See xclim.sdba.utils.extrapolate_qm() for details. Defaults to “constant”.

References

Cannon2015

Cannon, A. J., Sobie, S. R., & Murdock, T. Q. (2015). Bias correction of GCM precipitation by quantile mapping: How well do methods preserve changes in quantiles and extremes? Journal of Climate, 28(17), 6938–6959. https://doi.org/10.1175/JCLI-D-14-00754.1

class xclim.sdba.adjustment.EmpiricalQuantileMapping(*args, _trained=False, **kwargs)[source]

Bases: xclim.sdba.adjustment.TrainAdjust

Empirical Quantile Mapping bias-adjustment.

Adjustment factors are computed between the quantiles of ref and sim. Values of sim are matched to the corresponding quantiles of hist and corrected accordingly.

\[F^{-1}_{ref} (F_{hist}(sim))\]

where \(F\) is the cumulative distribution function (CDF) and mod stands for model data.

Parameters
  • Train step

  • nquantiles (int or 1d array of floats) – The number of quantiles to use. Two endpoints at 1e-6 and 1 - 1e-6 will be added. An array of quantiles [0, 1] can also be passed. Defaults to 20 quantiles (+2 endpoints).

  • kind ({‘+’, ‘’}*) – The adjustment kind, either additive or multiplicative. Defaults to “+”.

  • group (Union[str, Grouper]) – The grouping information. See xclim.sdba.base.Grouper for details. Default is “time”, meaning an single adjustment group along dimension “time”.

  • Adjust step

  • interp ({‘nearest’, ‘linear’, ‘cubic’}) – The interpolation method to use when interpolating the adjustment factors. Defaults to “nearset”.

  • extrapolation ({‘constant’, ‘nan’}) – The type of extrapolation to use. See xclim.sdba.utils.extrapolate_qm() for details. Defaults to “constant”.

References

Dequé, M. (2007). Frequency of precipitation and temperature extremes over France in an anthropogenic scenario: Model results and statistical correction according to observed values. Global and Planetary Change, 57(1–2), 16–26. https://doi.org/10.1016/j.gloplacha.2006.11.030

_allow_diff_calendars = False

class xclim.sdba.adjustment.ExtremeValues(*args, _trained=False, **kwargs)[source]

Bases: xclim.sdba.adjustment.TrainAdjust

Adjustment correction for extreme values.

The tail of the distribution of adjusted data is corrected according to the bias between the parametric Generalized Pareto distributions of the simulatated and reference data, [RRJF2021]. The distributions are composed of the maximal values of clusters of “large” values. With “large” values being those above cluster_thresh. Only extreme values, whose quantile within the pool of large values are above q_thresh, are re-adjusted. See Notes.

This adjustment method should be considered experimental and used with care.

Parameters
  • Train step

  • cluster_thresh (Quantity (str with units)) – The threshold value for defining clusters.

  • q_thresh (float) – The quantile of “extreme” values, [0, 1[. Defaults to 0.95.

  • ref_params (xr.DataArray, optional) – Distribution parameters to use instead of fitting a GenPareto distribution on ref.

  • Adjust step

  • scen (DataArray) – This is a second-order adjustment, so the adjust method needs the first-order adjusted timeseries in addition to the raw “sim”.

  • interp ({‘nearest’, ‘linear’, ‘cubic’}) – The interpolation method to use when interpolating the adjustment factors. Defaults to “linear”.

  • extrapolation ({‘constant’, ‘nan’}) – The type of extrapolation to use. See extrapolate_qm() for details. Defaults to “constant”.

  • frac (float) – Fraction where the cutoff happens between the original scen and the corrected one. See Notes, ]0, 1]. Defaults to 0.25.

  • power (float) – Shape of the correction strength, see Notes. Defaults to 1.0.

Notes

Extreme values are extracted from ref, hist and sim by finding all “clusters”, i.e. runs of consecutive values above cluster_thresh. The q_thresh`th percentile of these values is taken on `ref and hist and becomes thresh, the extreme value threshold. The maximal value of each cluster, if it exceeds that new threshold, is taken and Generalized Pareto distributions are fitted to them, for both ref and hist. The probabilities associated with each of these extremes in hist is used to find the corresponding value according to ref’s distribution. Adjustment factors are computed as the bias between those new extremes and the original ones.

In the adjust step, a Generalized Pareto distributions is fitted on the cluster-maximums of sim and it is used to associate a probability to each extreme, values over the thresh compute in the training, without the clustering. The adjustment factors are computed by interpolating the trained ones using these probabilities and the probabilities computed from hist.

Finally, the adjusted values (\(C_i\)) are mixed with the pre-adjusted ones (scen, \(D_i\)) using the following transition function:

\[V_i = C_i * \tau + D_i * (1 - \tau)\]

Where \(\tau\) is a function of sim’s extreme values (unadjusted, \(S_i\)) and of arguments frac (\(f\)) and power (\(p\)):

\[\tau = \left(\frac{1}{f}\frac{S - min(S)}{max(S) - min(S)}\right)^p\]

Code based on an internal Matlab source and partly ib the biascorrect_extremes function of the julia package [ClimateTools].

Because of limitations imposed by the lazy computing nature of the dask backend, it is not possible to know the number of cluster extremes in ref and hist at the moment the output data structure is created. This is why the code tries to estimate that number and usually overstimates it. In the training dataset, this translated into a quantile dimension that is too large and variables af and px_hist are assigned NaNs on extra elements. This has no incidence on the calculations themselves but requires more memory than is useful.

References

RRJF2021

Roy, P., Rondeau-Genesse, G., Jalbert, J., Fournier, É. 2021. Climate Scenarios of Extreme Precipitation Using a Combination of Parametric and Non-Parametric Bias Correction Methods. Submitted to Climate Services, April 2021.

ClimateTools

https://juliaclimate.github.io/ClimateTools.jl/stable/

class xclim.sdba.adjustment.LOCI(*args, _trained=False, **kwargs)[source]

Bases: xclim.sdba.adjustment.TrainAdjust

Local Intensity Scaling (LOCI) bias-adjustment.

This bias adjustment method is designed to correct daily precipitation time series by considering wet and dry days separately ([Schmidli2006]).

Multiplicative adjustment factors are computed such that the mean of hist matches the mean of ref for values above a threshold.

The threshold on the training target ref is first mapped to hist by finding the quantile in hist having the same exceedance probability as thresh in ref. The adjustment factor is then given by

\[s = \frac{\left \langle ref: ref \geq t_{ref} \right\rangle - t_{ref}}{\left \langle hist : hist \geq t_{hist} \right\rangle - t_{hist}}\]

In the case of precipitations, the adjustment factor is the ratio of wet-days intensity.

For an adjustment factor s, the bias-adjustment of sim is:

\[sim(t) = \max\left(t_{ref} + s \cdot (hist(t) - t_{hist}), 0\right)\]
Parameters
  • Train step

  • group (Union[str, Grouper]) – The grouping information. See xclim.sdba.base.Grouper for details. Default is “time”, meaning an single adjustment group along dimension “time”.

  • thresh (str) – The threshold in ref above which the values are scaled.

  • Adjust step

  • interp ({‘nearest’, ‘linear’, ‘cubic’}) – The interpolation method to use then interpolating the adjustment factors. Defaults to “linear”.

References

Schmidli2006

Schmidli, J., Frei, C., & Vidale, P. L. (2006). Downscaling from GCM precipitation: A benchmark for dynamical and statistical downscaling methods. International Journal of Climatology, 26(5), 679–689. DOI:10.1002/joc.1287

class xclim.sdba.adjustment.NpdfTransform(*args, _trained=False, **kwargs)[source]

Bases: xclim.sdba.adjustment.Adjust

N-dimensional probability density function transform.

This adjustment object combines both training and adjustent steps in the adjust class method.

A multivariate bias-adjustment algorithm described by [Cannon18], as part of the MBCn algorithm, based on a color-correction algorithm described by [Pitie05].

This algorithm in itself, when used with QuantileDeltaMapping, is NOT trend-preserving. The full MBCn algorithm includes a reordering step provided here by xclim.sdba.processing.reordering().

See notes for an explanation of the algorithm.

Parameters
  • base (BaseAdjustment) – An univariate bias-adjustment class. This is untested for anything else than QuantileDeltaMapping.

  • base_kws (dict, optional) – Arguments passed to the training of the univariate adjustment.

  • n_escore (int) – The number of elements to send to the escore function. The default, 0, means all elements are included. Pass -1 to skip computing the escore completely. Small numbers result in less significative scores, but the execution time goes up quickly with large values.

  • n_iter (int) – The number of iterations to perform. Defaults to 20.

  • pts_dim (str) – The name of the “multivariate” dimension. Defaults to “multivar”, which is the normal case when using xclim.sdba.base.stack_variables().

  • adj_kws (dict, optional) – Dictionary of arguments to pass to the adjust method of the univariate adjustment.

  • rot_matrices (xr.DataArray, optional) – The rotation matrices as a 3D array (‘iterations’, <pts_dim>, <anything>), with shape (n_iter, <N>, <N>). If left empty, random rotation matrices will be automatically generated.

Notes

The historical reference (\(T\), for “target”), simulated historical (\(H\)) and simulated projected (\(S\)) datasets are constructed by stacking the timeseries of N variables together. The algoriths goes into the following steps:

  1. Rotate the datasets in the N-dimensional variable space with \(\mathbf{R}\), a random rotation NxN matrix.

..math

\tilde{\mathbf{T}} = \mathbf{T}\mathbf{R} \\
\tilde{\mathbf{H}} = \mathbf{H}\mathbf{R} \\
\tilde{\mathbf{S}} = \mathbf{S}\mathbf{R}

2. An univariate bias-adjustment \(\mathcal{F}\) is used on the rotated datasets. The adjustments are made in additive mode, for each variable \(i\).

\[\hat{\mathbf{H}}_i, \hat{\mathbf{S}}_i = \mathcal{F}\left(\tilde{\mathbf{T}}_i, \tilde{\mathbf{H}}_i, \tilde{\mathbf{S}}_i\right)\]
  1. The bias-adjusted datasets are rotated back.

\[\begin{split}\mathbf{H}' = \hat{\mathbf{H}}\mathbf{R} \\ \mathbf{S}' = \hat{\mathbf{S}}\mathbf{R}\end{split}\]

These three steps are repeated a certain number of times, prescribed by argument n_iter. At each iteration, a new random rotation matrix is generated.

The original algorithm ([Pitie05]), stops the iteration when some distance score converges. Following [Cannon18] and the MBCn implementation in [CannonR], we instead fix the number of iterations.

As done by [Cannon18], the distance score chosen is the “Energy distance” from [SkezelyRizzo] (see xclim.sdba.processing.escore()).

The random matrices are generated following a method laid out by [Mezzadri].

This is only part of the full MBCn algorithm, see Fourth example : Multivariate bias-adjustment with multiple steps - Cannon 2018 for an example on how to replicate the full method with xclim. This includes a standardization of the simulated data beforehand, an initial univariate adjustment and the reordering of those adjusted series according to the rank structure of the output of this algorithm.

References

Cannon18

Cannon, A. J. (2018). Multivariate quantile mapping bias correction: An N-dimensional probability density function transform for climate model simulations of multiple variables. Climate Dynamics, 50(1), 31–49. https://doi.org/10.1007/s00382-017-3580-6

CannonR

https://CRAN.R-project.org/package=MBC

Mezzadri

Mezzadri, F. (2006). How to generate random matrices from the classical compact groups. arXiv preprint math-ph/0609050.

Pitie05(1,2)

Pitie, F., Kokaram, A. C., & Dahyot, R. (2005). N-dimensional probability density function transfer and its application to color transfer. Tenth IEEE International Conference on Computer Vision (ICCV’05) Volume 1, 2, 1434-1439 Vol. 2. https://doi.org/10.1109/ICCV.2005.166

SkezelyRizzo

Szekely, G. J. and Rizzo, M. L. (2004) Testing for Equal Distributions in High Dimension, InterStat, November (5)

class xclim.sdba.adjustment.PrincipalComponents(*args, _trained=False, **kwargs)[source]

Bases: xclim.sdba.adjustment.TrainAdjust

Principal component adjustment.

This bias-correction method maps model simulation values to the observation space through principal components ([hnilica2017]). Values in the simulation space (multiple variables, or multiple sites) can be thought of as coordinates along axes, such as variable, temperature, etc. Principal components (PC) are a linear combinations of the original variables where the coefficients are the eigenvectors of the covariance matrix. Values can then be expressed as coordinates along the PC axes. The method makes the assumption that bias-corrected values have the same coordinates along the PC axes of the observations. By converting from the observation PC space to the original space, we get bias corrected values. See notes for a mathematical explanation.

Note that principal components is meant here as the algebraic operation defining a coordinate system based on the eigenvectors, not statistical principal component analysis.

Parameters
  • group (Union[str, Grouper]) – The main dimension and grouping information. See Notes. See xclim.sdba.base.Grouper for details. The adjustment will be performed on each group independently. Default is “time”, meaning an single adjustment group along dimension “time”.

  • best_orientation ({‘simple’, ‘full’}) – Which method to use when searching for the best principal component orientation. See best_pc_orientation_simple() and best_pc_orientiation_full(). “full” is more precise, but it is much slower.

  • crd_dim (str) – The data dimension along which the multiple simulation space dimensions are taken. For a multivariate adjustment, this usually is “multivar”, as returned by sdba.stack_variables. For a multisite adjustment, this should be the spatial dimension. The training algorithm currently doesn’t support any chunking along either crd_dim. group.dim and group.add_dims.

Notes

The input data is understood as a set of N points in a \(M\)-dimensional space.

  • \(M\) is taken along crd_dim.

  • \(N\) is taken along the dimensions given through group : (the main dim but also, if requested, the add_dims and window).

The principal components (PC) of hist and ref are used to defined new coordinate systems, centered on their respective means. The training step creates a matrix defining the transformation from hist to ref:

\[scen = e_{R} + \mathrm{\mathbf{T}}(sim - e_{H})\]

Where:

\[\mathrm{\mathbf{T}} = \mathrm{\mathbf{R}}\mathrm{\mathbf{H}}^{-1}\]

\(\mathrm{\mathbf{R}}\) is the matrix transforming from the PC coordinates computed on ref to the data coordinates. Similarly, \(\mathrm{\mathbf{H}}\) is transform from the hist PC to the data coordinates (\(\mathrm{\mathbf{H}}\) is the inverse transformation). \(e_R\) and \(e_H\) are the centroids of the ref and hist distributions respectively. Upon running the adjust step, one may decide to use \(e_S\), the centroid of the sim distribution, instead of \(e_H\).

References

hnilica2017

Hnilica, J., Hanel, M. and Pš, V. (2017), Multisite bias correction of precipitation data from regional climate models. Int. J. Climatol., 37: 2934-2946. https://doi.org/10.1002/joc.4890

class xclim.sdba.adjustment.QuantileDeltaMapping(*args, _trained=False, **kwargs)[source]

Bases: xclim.sdba.adjustment.EmpiricalQuantileMapping

Quantile Delta Mapping bias-adjustment.

Adjustment factors are computed between the quantiles of ref and hist. Quantiles of sim are matched to the corresponding quantiles of hist and corrected accordingly.

\[sim\frac{F^{-1}_{ref}\left[F_{sim}(sim)\right]}{F^{-1}_{hist}\left[F_{sim}(sim)\right]}\]

where \(F\) is the cumulative distribution function (CDF). This equation is valid for multiplicative adjustment. The algorithm is based on the “QDM” method of [Cannon2015].

Parameters
  • Train step

  • nquantiles (int or 1d array of floats) – The number of quantiles to use. Two endpoints at 1e-6 and 1 - 1e-6 will be added. An array of quantiles [0, 1] can also be passed. Defaults to 20 quantiles (+2 endpoints).

  • kind ({‘+’, ‘’}*) – The adjustment kind, either additive or multiplicative. Defaults to “+”.

  • group (Union[str, Grouper]) – The grouping information. See xclim.sdba.base.Grouper for details. Default is “time”, meaning an single adjustment group along dimension “time”.

  • Adjust step

  • interp ({‘nearest’, ‘linear’, ‘cubic’}) – The interpolation method to use when interpolating the adjustment factors. Defaults to “nearest”.

  • extrapolation ({‘constant’, ‘nan’}) – The type of extrapolation to use. See xclim.sdba.utils.extrapolate_qm() for details. Defaults to “constant”.

  • Extra diagnostics

  • —————–

  • In adjustment

  • quantiles (The quantile of each value of sim. The adjustment factor is interpolated using this as the “quantile” axis on ds.af.)

References

Cannon2015

Cannon, A. J., Sobie, S. R., & Murdock, T. Q. (2015). Bias correction of GCM precipitation by quantile mapping: How well do methods preserve changes in quantiles and extremes? Journal of Climate, 28(17), 6938–6959. https://doi.org/10.1175/JCLI-D-14-00754.1

class xclim.sdba.adjustment.Scaling(*args, _trained=False, **kwargs)[source]

Bases: xclim.sdba.adjustment.TrainAdjust

Scaling bias-adjustment.

Simple bias-adjustment method scaling variables by an additive or multiplicative factor so that the mean of hist matches the mean of ref.

Parameters
  • Train step

  • group (Union[str, Grouper]) – The grouping information. See xclim.sdba.base.Grouper for details. Default is “time”, meaning an single adjustment group along dimension “time”.

  • kind ({‘+’, ‘’}*) – The adjustment kind, either additive or multiplicative. Defaults to “+”.

  • Adjust step

  • interp ({‘nearest’, ‘linear’, ‘cubic’}) – The interpolation method to use then interpolating the adjustment factors. Defaults to “nearest”.

Pre and post processing

xclim.sdba.processing.adapt_freq(ref: xarray.DataArray, sim: xarray.DataArray, *, group: Union[xclim.sdba.base.Grouper, str], thresh: str = '0 mm d-1') xarray.Dataset[source]

Adapt frequency of values under thresh of sim, in order to match ref.

This is useful when the dry-day frequency in the simulations is higher than in the references. This function will create new non-null values for sim/hist, so that adjustment factors are less wet-biased. Based on [Themessl2012].

Parameters
  • ds (xr.Dataset) – With variables : “ref”, Target/reference data, usually observed data. and “sim”, Simulated data.

  • dim (str) – Dimension name.

  • group (Union[str, Grouper]) – Grouping information, see base.Grouper

  • thresh (str) – Threshold below which values are considered zero, a quantity with units.

Returns

  • sim_adj (xr.DataArray) – Simulated data with the same frequency of values under threshold than ref. Adjustment is made group-wise.

  • pth (xr.DataArray) – For each group, the smallest value of sim that was not frequency-adjusted. All values smaller were either left as zero values or given a random value between thresh and pth. NaN where frequency adaptation wasn’t needed.

  • dP0 (xr.DataArray) – For each group, the percentage of values that were corrected in sim.

Notes

With \(P_0^r\) the frequency of values under threshold \(T_0\) in the reference (ref) and \(P_0^s\) the same for the simulated values, \(\\Delta P_0 = \\frac{P_0^s - P_0^r}{P_0^s}\), when positive, represents the proportion of values under \(T_0\) that need to be corrected.

The correction replaces a proportion \(\\Delta P_0\) of the values under \(T_0\) in sim by a uniform random number between \(T_0\) and \(P_{th}\), where \(P_{th} = F_{ref}^{-1}( F_{sim}( T_0 ) )\) and F(x) is the empirical cumulative distribution function (CDF).

References

Themessl2012

Themeßl et al. (2012), Empirical-statistical downscaling and error correction of regional climate models and its impact on the climate change signal, Climatic Change, DOI 10.1007/s10584-011-0224-4.

xclim.sdba.processing.construct_moving_yearly_window(da: xarray.Dataset, window: int = 21, step: int = 1, dim: str = 'movingwin')[source]

Construct a moving window DataArray.

Stacks windows of da in a new ‘movingwin’ dimension. Windows are always made of full years, so calendar with non-uniform year lengths are not supported.

Windows are constructed starting at the beginning of da, if number of given years is not a multiple of step, then the last year(s) will be missing as a supplementary window would be incomplete.

Parameters
  • da (xr.Dataset) – A DataArray with a time dimension.

  • window (int) – The length of the moving window as a number of years.

  • step (int) – The step between each window as a number of years.

  • dim (str) – The new dimension name. If given, must also be given to unpack_moving_yearly_window.

Returns

xr.DataArray – A DataArray with a new movingwin dimension and a time dimension with a length of 1 window. This assumes downstream algorithms do not make use of the _absolute_ year of the data. The correct timeseries can be reconstructed with unpack_moving_yearly_window(). The coordinates of movingwin are the first date of the windows.

xclim.sdba.processing.escore(tgt: xarray.DataArray, sim: xarray.DataArray, dims: Sequence[str] = ('variables', 'time'), N: int = 0, scale: bool = False) xarray.DataArray[source]

Energy score, or energy dissimilarity metric, based on [SkezelyRizzo] and [Cannon18].

Parameters
  • tgt (xr.DataArray) – Target observations.

  • sim (xr.DataArray) – Candidate observations. Must have the same dimensions as tgt.

  • dims (sequence of 2 strings) – The name of the dimensions along which the variables and observation points are listed. tgt and sim can have different length along the second one, but must be equal along the first one. The result will keep all other dimensions.

  • N (int) – If larger than 0, the number of observations to use in the score computation. The points are taken evenly distributed along obs_dim.

  • scale (bool) – Whether to scale the data before computing the score. If True, both arrays as scaled according to the mean and standard deviation of tgt along obs_dim. (std computed with ddof=1 and both statistics excluding NaN values).

Returns

xr.DataArray – e-score with dimensions not in dims.

Notes

Explanation adapted from the “energy” R package documentation. The e-distance between two clusters \(C_i\), \(C_j\) (tgt and sim) of size \(n_i,n_j\) proposed by Skezely and Rizzo (2004) is defined by:

\[e(C_i,C_j) = \frac{1}{2}\frac{n_i n_j}{n_i + n_j} \left[2 M_{ij} − M_{ii} − M_{jj}\right]\]

where

\[M_{ij} = \frac{1}{n_i n_j} \sum_{p = 1}^{n_i} \sum{q = 1}^{n_j} \left\Vert X_{ip} − X{jq} \right\Vert.\]

\(\Vert\cdot\Vert\) denotes Euclidean norm, \(X_{ip}\) denotes the p-th observation in the i-th cluster.

The input scaling and the factor \(\frac{1}{2}\) in the first equation are additions of [Cannon18] to the metric. With that factor, the test becomes identical to the one defined by [BaringhausFranz].

References

BaringhausFranz

Baringhaus, L. and Franz, C. (2004) On a new multivariate two-sample test, Journal of Multivariate Analysis, 88(1), 190–206. https://doi.org/10.1016/s0047-259x(03)00079-4

Cannon18

Cannon, A. J. (2018). Multivariate quantile mapping bias correction: An N-dimensional probability density function transform for climate model simulations of multiple variables. Climate Dynamics, 50(1), 31–49. https://doi.org/10.1007/s00382-017-3580-6

SkezelyRizzo

Skezely, G. J. and Rizzo, M. L. (2004) Testing for Equal Distributions in High Dimension, InterStat, November (5)

xclim.sdba.processing.from_additive_space(data: xarray.DataArray, lower_bound: Optional[str] = None, upper_bound: Optional[str] = None, trans: Optional[str] = None, units: Optional[str] = None)[source]

Transform back to the physical space a variable that was transformed with to_addtitive_space.

Based on [AlavoineGrenier]. If parameters are not present on the attributes of the data, they must be all given are arguments.

Parameters
  • data (xr.DataArray) – A variable that was transform by to_additive_space().

  • lower_bound (str, optional) – The smallest physical value of the variable, as a Quantity string. The final data will have no value smaller or equal to this bound. If None (default), the sdba_transform_lower attribute is looked up on data.

  • upper_bound (str, optional) – The largest physical value of the variable, as a Quantity string. Only relevant for the logit transformation. The final data will have no value larger or equal to this bound. If None (default), the sdba_transform_upper attribute is looked up on data.

  • trans ({‘log’, ‘logit’}, optional) – The transformation to use. See notes. If None (the default), the sdba_transform attribute is looked up on data.

  • units (str, optional) – The units of the data before transformation to the additive space. If None (the default), the sdba_transform_units attribute is looked up on data.

Returns

xr.DataArray – The physical variable. Attributes are conserved, even if some might be incorrect. Except units which are taken from sdba_transform_units if available. All sdba_transform* attributes are deleted.

Notes

Given a variable that is not usable in an additive adjustment, to_additive_space() applied a transformation to a space where additive methods are sensible. Given \(Y\) the transformed variable, \(b_-\) the lower physical bound of that variable and \(b_+\) the upper physical bound, two back-transformations are currently implemented to get \(X\), the physical variable.

  • log

    \[X = e^{Y} + b_-\]
  • logit

    \[X' = \frac{1}{1 + e^{-Y}} X = X * (b_+ - b_-) + b_-\]

See also

to_additive_space

for the original transformation.

References

AlavoineGrenier

Alavoine M., and Grenier P. (under review) The distinct problems of physical inconsistency and of multivariate bias potentially involved in the statistical adjustment of climate simulations. International Journal of Climatology, Manuscript ID: JOC-21-0789, submitted on September 19th 2021. (Preprint https://doi.org/10.31223/X5C34C)

xclim.sdba.processing.jitter(x: xarray.DataArray, lower: Optional[str] = None, upper: Optional[str] = None, minimum: Optional[str] = None, maximum: Optional[str] = None) xarray.DataArray[source]

Replaces values under a threshold and values above another by a uniform random noise.

Do not confuse with R’s jitter, which adds uniform noise instead of replacing values.

Parameters
  • x (xr.DataArray) – Values.

  • lower (str) – Threshold under which to add uniform random noise to values, a quantity with units. If None, no jittering is performed on the lower end.

  • upper (str) – Threshold over which to add uniform random noise to values, a quantity with units. If None, no jittering is performed on the upper end.

  • minimum (str) – Lower limit (excluded) for the lower end random noise, a quantity with units. If None but lower is not None, 0 is used.

  • maximum (str) – Upper limit (excluded) for the upper end random noise, a quantity with units. If upper is not None, it must be given.

Returns

xr.DataArray – Same as x but values < lower are replaced by a uniform noise in range (minimum, lower) and values >= upper are replaced by a uniform noise in range [upper, maximum). The two noise distributions are independent.

xclim.sdba.processing.jitter_over_thresh(x: xarray.DataArray, thresh: str, upper_bnd: str) xarray.DataArray[source]

Replace values greater than threshold by a uniform random noise.

Do not confuse with R’s jitter, which adds uniform noise instead of replacing values.

Parameters
  • x (xr.DataArray) – Values.

  • thresh (str) – Threshold over which to add uniform random noise to values, a quantity with units.

  • upper_bnd (str) – Maximum possible value for the random noise, a quantity with units.

Returns

xr.DataArray

Notes

If thresh is low, this will change the mean value of x.

xclim.sdba.processing.jitter_under_thresh(x: xarray.DataArray, thresh: str) xarray.DataArray[source]

Replace values smaller than threshold by a uniform random noise.

Do not confuse with R’s jitter, which adds uniform noise instead of replacing values.

Parameters
  • x (xr.DataArray) – Values.

  • thresh (str) – Threshold under which to add uniform random noise to values, a quantity with units.

Returns

xr.DataArray

Notes

If thresh is high, this will change the mean value of x.

xclim.sdba.processing.normalize(data: xarray.DataArray, norm: Optional[xarray.DataArray] = None, *, group: Union[xclim.sdba.base.Grouper, str], kind: str = '+') xarray.Dataset[source]

Normalize an array by removing its mean.

Normalization if performed group-wise and according to kind.

Parameters
  • data (xr.DataArray) – The variable to normalize.

  • norm (xr.DataArray, optional) – If present, it is used instead of computing the norm again.

  • group (Union[str, Grouper]) – Grouping information. See xclim.sdba.base.Grouper for details..

  • kind ({‘+’, ‘’}*) – If kind is “+”, the mean is subtracted from the mean and if it is ‘*’, it is divided from the data.

Returns

  • xr.DataArray – Groupwise anomaly.

  • norm (xr.DataArray) – Mean over each group.

xclim.sdba.processing.reordering(ref: xarray.DataArray, sim: xarray.DataArray, group: str = 'time') xarray.Dataset[source]

Reorders data in sim following the order of ref.

The rank structure of ref is used to reorder the elements of sim along dimension “time”, optionally doing the operation group-wise.

Parameters
  • sim (xr.DataArray) – Array to reorder.

  • ref (xr.DataArray) – Array whose rank order sim should replicate.

  • group (str) – Grouping information. See xclim.sdba.base.Grouper for details.

Returns

  • xr.Dataset – sim reordered according to ref’s rank order.

  • Reference

  • ———

  • .. [Cannon18] Cannon, A. J. (2018). Multivariate quantile mapping bias correction (An N-dimensional probability density function transform for climate model simulations of multiple variables. Climate Dynamics, 50(1), 31–49. https://doi.org/10.1007/s00382-017-3580-6)

xclim.sdba.processing.stack_variables(ds: xarray.Dataset, rechunk: bool = True, dim: str = 'multivar')[source]

Stack different variables of a dataset into a single DataArray with a new “variables” dimension.

Variable attributes are all added as lists of attributes to the new coordinate, prefixed with “_”. Variables are concatenated in the new dimension in alphabetical order, to ensure coherent behaviour with different datasets.

Parameters
  • ds (xr.Dataset) – Input dataset.

  • rechunk (bool) – If True (default), dask arrays are rechunked with variables : -1.

  • dim (str) – Name of dimension along which variables are indexed.

Returns

xr.DataArray – The transformed variable. Attributes are conserved, even if some might be incorrect. Except units, which are replaced with “”. Old units are stored in sdba_transformation_units. A sdba_transform attribute is added, set to the transformation method. sdba_transform_lower and sdba_transform_upper are also set if the requested bounds are different from the defaults.

Array with variables stacked along dim dimension. Units are set to “”.

xclim.sdba.processing.standardize(da: xarray.DataArray, mean: Optional[xarray.DataArray] = None, std: Optional[xarray.DataArray] = None, dim: str = 'time') Tuple[Union[xarray.DataArray, xarray.Dataset], xarray.DataArray, xarray.DataArray][source]

Standardize a DataArray by centering its mean and scaling it by its standard deviation.

Either of both of mean and std can be provided if need be.

Returns the standardized data, the mean and the standard deviation.

xclim.sdba.processing.to_additive_space(data: xarray.DataArray, lower_bound: str, upper_bound: Optional[str] = None, trans: str = 'log')[source]

Transform a non-additive variable into an additive space by the means of a log or logit transformation.

Based on [AlavoineGrenier].

Parameters
  • data (xr.DataArray) – A variable that can’t usually be bias-adusted by additive methods.

  • lower_bound (str) – The smallest physical value of the variable, excluded, as a Quantity string. The data should only have values strictly larger than this bound.

  • upper_bound (str, optional) – The largest physical value of the variable, excluded, as a Quantity string. Only relevant for the logit transformation. The data should only have values strictly smaller than this bound.

  • trans ({‘log’, ‘logit’}) – The transformation to use. See notes.

Notes

Given a variable that is not usable in an additive adjustment, this apply a transformation to a space where additive methods are sensible. Given \(X\) the variable, \(b_-\) the lower physical bound of that variable and \(b_+\) the upper physical bound, two transformations are currently implemented to get \(Y\), the additive-ready variable. \(\ln\) is the natural logarithm.

  • log

    \[Y = \ln\left( X - b_- \right)\]

    Usually used for variables with only a lower bound, like precipitation (pr, prsn, etc) and daily temperature range (dtr). Both have a lower bound of 0.

  • logit

    \[X' = (X - b_-) / (b_+ - b_-) Y = \ln\left(\frac{X'}{1 - X'} \right)\]

    Usually used for variables with both a lower and a upper bound, like relative and specific humidity, cloud cover fraction, etc.

This will thus produce Infinity and NaN values where \(X == b_-\) or \(X == b_+\). We recommend using jitter_under_thresh() and jitter_over_thresh() to remove those issues.

See also

from_additive_space

for the inverse transformation.

jitter_under_thresh

Remove values exactly equal to the lower bound.

jitter_over_thresh

Remove values exactly equal to the upper bound.

References

AlavoineGrenier

Alavoine M., and Grenier P. (under review) The distinct problems of physical inconsistency and of multivariate bias potentially involved in the statistical adjustment of climate simulations. International Journal of Climatology, Manuscript ID: JOC-21-0789, submitted on September 19th 2021. (Preprint https://doi.org/10.31223/X5C34C)

xclim.sdba.processing.uniform_noise_like(da: xarray.DataArray, low: float = 1e-06, high: float = 0.001) xarray.DataArray[source]

Return a uniform noise array of the same shape as da.

Noise is uniformly distributed between low and high. Alternative method to jitter_under_thresh for avoiding zeroes.

xclim.sdba.processing.unpack_moving_yearly_window(da: xarray.DataArray, dim: str = 'movingwin')[source]

Unpack a constructed moving window dataset to a normal timeseries, only keeping the central data.

Unpack DataArrays created with construct_moving_yearly_window() and recreate a timeseries data. Only keeps the central non-overlapping years. The final timeseries will be (window - step) years shorter than the initial one.

The window length and window step are inferred from the coordinates.

Parameters
xclim.sdba.processing.unstack_variables(da: xarray.DataArray, dim: Optional[str] = None)[source]

Unstack a DataArray created by stack_variables to a dataset.

Parameters
  • da (xr.DataArray) – Array holding different variables along dim dimension.

  • dim (str) – Name of dimension along which the variables are stacked. If not specified (default), dim is inferred from attributes of the coordinate.

Returns

xr.Dataset – Dataset holding each variable in an individual DataArray.

xclim.sdba.processing.unstandardize(da: xarray.DataArray, mean: xarray.DataArray, std: xarray.DataArray)[source]

Rescale a standardized array by performing the inverse operation of standardize.

Detrending objects

class xclim.sdba.detrending.LoessDetrend(group='time', kind='+', f=0.2, niter=1, d=0, weights='tricube', equal_spacing=None)[source]

Bases: xclim.sdba.detrending.BaseDetrend

Detrend time series using a LOESS regression.

The fit is a piecewise linear regression. For each point, the contribution of all neighbors is weighted by a bell-shaped curve (gaussian) with parameters sigma (std). The x-coordinate of the DataArray is scaled to [0,1] before the regression is computed.

Parameters
  • group (Union[str, Grouper]) – The grouping information. See xclim.sdba.base.Grouper for details. The fit is performed along the group’s main dim.

  • kind ({’’, ‘+’}*) – The way the trend is removed or added, either additive or multiplicative.

  • d ([0, 1]) – Order of the local regression. Only 0 and 1 currently implemented.

  • f (float) – Parameter controlling the span of the weights, between 0 and 1.

  • niter (int) – Number of robustness iterations to execute.

  • weights ([“tricube”, “gaussian”]) – Shape of the weighting function: “tricube” : a smooth top-hat like curve, f gives the span of non-zero values. “gaussian” : a gaussian curve, f gives the span for 95% of the values.

Notes

LOESS smoothing is computationally expensive. As it relies on a loop on gridpoints, it can be useful to use smaller than usual chunks. Moreover, it suffers from heavy boundary effects. As a rule of thumb, the outermost N * f/2 points should be considered dubious. (N is the number of points along each group)

class xclim.sdba.detrending.MeanDetrend(*, group: Union[xclim.sdba.base.Grouper, str] = 'time', kind: str = '+', **kwargs)[source]

Bases: xclim.sdba.detrending.BaseDetrend

Simple detrending removing only the mean from the data, quite similar to normalizing.

class xclim.sdba.detrending.NoDetrend(*, group: Union[xclim.sdba.base.Grouper, str] = 'time', kind: str = '+', **kwargs)[source]

Bases: xclim.sdba.detrending.BaseDetrend

Convenience class for polymorphism. Does nothing.

class xclim.sdba.detrending.PolyDetrend(group='time', kind='+', degree=4, preserve_mean=False)[source]

Bases: xclim.sdba.detrending.BaseDetrend

Detrend time series using a polynomial regression.

Parameters
  • group (Union[str, Grouper]) – The grouping information. See xclim.sdba.base.Grouper for details. The fit is performed along the group’s main dim.

  • kind ({’’, ‘+’}*) – The way the trend is removed or added, either additive or multiplicative.

  • degree (int) – The order of the polynomial to fit.

  • preserve_mean (bool) – Whether to preserve the mean when de/re-trending. If True, the trend has its mean removed before it is used.

class xclim.sdba.detrending.RollingMeanDetrend(group='time', kind='+', win=30, weights=None, min_periods=None)[source]

Bases: xclim.sdba.detrending.BaseDetrend

Detrend time series using a rolling mean.

Parameters
  • group (Union[str, Grouper]) – The grouping information. See xclim.sdba.base.Grouper for details. The fit is performed along the group’s main dim.

  • kind ({’’, ‘+’}*) – The way the trend is removed or added, either additive or multiplicative.

  • win (int) – The size of the rolling window. Units are the steps of the grouped data, which means this detrending is best use with either group=’time’ or group=’time.dayofyear’. Other grouping will have large jumps included within the windows and :py`:class:LoessDetrend might offer a better solution.

  • weights (sequence of floats, optional) – Sequence of length win. Defaults to None, which means a flat window.

  • min_periods (int, optional) – Minimum number of observations in window required to have a value, otherwise the result is NaN. See xarray.DataArray.rolling(). Defaults to None, which sets it equal to win. Setting both weights and this is not implemented yet.

Notes

As for the LoessDetrend detrending, important boundary effects are to be expected.

sdba utilities

xclim.sdba.utils.add_cyclic_bounds(da: xarray.DataArray, att: str, cyclic_coords: bool = True) Union[xarray.DataArray, xarray.Dataset][source]

Reindex an array to include the last slice at the beginning and the first at the end.

This is done to allow interpolation near the end-points.

Parameters
  • da (Union[xr.DataArray, xr.Dataset]) – An array

  • att (str) – The name of the coordinate to make cyclic

  • cyclic_coords (bool) – If True, the coordinates are made cyclic as well, if False, the new values are guessed using the same step as their neighbour.

Returns

Union[xr.DataArray, xr.Dataset] – da but with the last element along att prepended and the last one appended.

xclim.sdba.utils.apply_correction(x: xarray.DataArray, factor: xarray.DataArray, kind: Optional[str] = None) xarray.DataArray[source]

Apply the additive or multiplicative correction/adjustment factors.

If kind is not given, default to the one stored in the “kind” attribute of factor.

xclim.sdba.utils.best_pc_orientation_full(R: numpy.ndarray, Hinv: numpy.ndarray, Rmean: numpy.ndarray, Hmean: numpy.ndarray, hist: numpy.ndarray) numpy.ndarray[source]

Return best orientation vector for A according to the method of Alavoine et al. (2021, preprint).

Eigenvectors returned by pc_matrix do not have a defined orientation. Given an inverse transform Hinv, a transform R, the actual and target origins Hmean and Rmean and the matrix of training observations hist, this computes a scenario for all possible orientations and return the orientation that maximizes the Spearman correlation coefficient of all variables. The correlation is computed for each variable individually, then averaged.

This trick is explained in [alavoine2021]. See documentation of sdba.adjustment.PrincipalComponentAdjustment().

Parameters
  • R (np.ndarray) – MxM Matrix defining the final transformation.

  • Hinv (np.ndarray) – MxM Matrix defining the (inverse) first transformation.

  • Rmean (np.ndarray) – M vector defining the target distribution center point.

  • Hmean (np.ndarray) – M vector defining the original distribution center point.

  • hist (np.ndarray) – MxN matrix of all training observations of the M variables/sites.

Returns

np.ndarray – M vector of orientation correction (1 or -1).

References

alavoine2021

Alavoine, M., & Grenier, P. (2021). The distinct problems of physical inconsistency and of multivariate bias potentially involved in the statistical adjustment of climate simulations. https://eartharxiv.org/repository/view/2876/

xclim.sdba.utils.best_pc_orientation_simple(R: numpy.ndarray, Hinv: numpy.ndarray, val: float = 1000) numpy.ndarray[source]

Return best orientation vector according to a simple test.

Eigenvectors returned by pc_matrix do not have a defined orientation. Given an inverse transform Hinv and a transform R, this returns the orientation minimizing the projected distance for a test point far from the origin.

This trick is inspired by the one exposed in [hnilica2017]. For each possible orientation vector, the test point is reprojected and the distance from the original point is computed. The orientation minimizing that distance is chosen. See documentation of sdba.adjustment.PrincipalComponentAdjustment.

Parameters
  • R (np.ndarray) – MxM Matrix defining the final transformation.

  • Hinv (np.ndarray) – MxM Matrix defining the (inverse) first transformation.

  • val (float) – The coordinate of the test point (same for all axes). It should be much greater than the largest furthest point in the array used to define B.

Returns

np.ndarray – Mx1 vector of orientation correction (1 or -1).

References

hnilica2017

Hnilica, J., Hanel, M. and Pš, V. (2017), Multisite bias correction of precipitation data from regional climate models. Int. J. Climatol., 37: 2934-2946. https://doi.org/10.1002/joc.4890

xclim.sdba.utils.broadcast(grouped: xarray.DataArray, x: xarray.DataArray, *, group: Union[str, xclim.sdba.base.Grouper] = 'time', interp: str = 'nearest', sel: Optional[Mapping[str, xarray.DataArray]] = None) xarray.DataArray[source]

Broadcast a grouped array back to the same shape as a given array.

Parameters
  • grouped (xr.DataArray) – The grouped array to broadcast like x.

  • x (xr.DataArray) – The array to broadcast grouped to.

  • group (Union[str, Grouper]) – Grouping information. See xclim.sdba.base.Grouper for details.

  • interp ({‘nearest’, ‘linear’, ‘cubic’}) – The interpolation method to use,

  • sel (Mapping[str, xr.DataArray]) – Mapping of grouped coordinates to x coordinates (other than the grouping one).

Returns

xr.DataArray

xclim.sdba.utils.copy_all_attrs(ds: Union[xarray.Dataset, xarray.DataArray], ref: Union[xarray.Dataset, xarray.DataArray])[source]

Copies all attributes of ds to ref, including attributes of shared coordinates, and variables in the case of Datasets.

xclim.sdba.utils.ecdf(x: xarray.DataArray, value: float, dim: str = 'time') xarray.DataArray[source]

Return the empirical CDF of a sample at a given value.

Parameters
  • x (array) – Sample.

  • value (float) – The value within the support of x for which to compute the CDF value.

  • dim (str) – Dimension name.

Returns

xr.DataArray – Empirical CDF.

xclim.sdba.utils.ensure_longest_doy(func: Callable) Callable[source]

Ensure that selected day is the longest day of year for x and y dims.

xclim.sdba.utils.equally_spaced_nodes(n: int, eps: Optional[float] = 0.0001) numpy.array[source]

Return nodes with n equally spaced points within [0, 1] plus two end-points.

Parameters
  • n (int) – Number of equally spaced nodes.

  • eps (float, None) – Distance from 0 and 1 of end nodes. If None, do not add endpoints.

Returns

np.array – Nodes between 0 and 1.

Notes

For n=4, eps=0 : 0—x——x——x——x—1

xclim.sdba.utils.get_clusters(data: xarray.DataArray, u1, u2, dim: str = 'time') xarray.Dataset[source]

Get cluster count, maximum and position along a given dim.

See get_clusters_1d. Used by adjustment.ExtremeValues.

Parameters
  • data (1D ndarray) – Values to get clusters from.

  • u1 (float) – Extreme value threshold, at least one value in the cluster must exceed this.

  • u2 (float) – Cluster threshold, values above this can be part of a cluster.

  • dim (str) – Dimension name.

Returns

xr.Dataset

With variables,
  • nclusters : Number of clusters for each point (with dim reduced), int

  • start : First index in the cluster (dim reduced, new cluster), int

  • end : Last index in the cluster, inclusive (dim reduced, new cluster), int

  • maxpos : Index of the maximal value within the cluster (dim reduced, new cluster), int

  • maximum : Maximal value within the cluster (dim reduced, new cluster), same dtype as data.

For start, end and maxpos, -1 means NaN and should always correspond to a NaN in maximum. The length along cluster is half the size of “dim”, the maximal theoretical number of clusters.

xclim.sdba.utils.get_clusters_1d(data: numpy.ndarray, u1: float, u2: float) Tuple[numpy.array, numpy.array, numpy.array, numpy.array][source]

Get clusters of a 1D array.

A cluster is defined as a sequence of values larger than u2 with at least one value larger than u1.

Parameters
  • data (1D ndarray) – Values to get clusters from.

  • u1 (float) – Extreme value threshold, at least one value in the cluster must exceed this.

  • u2 (float) – Cluster threshold, values above this can be part of a cluster.

Returns

(np.array, np.array, np.array, np.array)

References

getcluster of Extremes.jl (read on 2021-04-20) https://github.com/jojal5/Extremes.jl

xclim.sdba.utils.get_correction(x: xarray.DataArray, y: xarray.DataArray, kind: str) xarray.DataArray[source]

Return the additive or multiplicative correction/adjustment factors.

xclim.sdba.utils.interp_on_quantiles(newx: xarray.DataArray, xq: xarray.DataArray, yq: xarray.DataArray, *, group: Union[str, xclim.sdba.base.Grouper] = 'time', method: str = 'linear', extrapolation: str = 'constant')[source]

Interpolate values of yq on new values of x.

Interpolate in 2D with griddata() if grouping is used, in 1D otherwise, with interp1d. Any NaNs in xq or yq are removed from the input map. Similarly, NaNs in newx are left NaNs.

Parameters
  • newx (xr.DataArray) – The values at which to evaluate yq. If group has group information, new should have a coordinate with the same name as the group name In that case, 2D interpolation is used.

  • xq, yq (xr.DataArray) – Coordinates and values on which to interpolate. The interpolation is done along the “quantiles” dimension if group has no group information. If it does, interpolation is done in 2D on “quantiles” and on the group dimension.

  • group (Union[str, Grouper]) – The dimension and grouping information. (ex: “time” or “time.month”). Defaults to “time”.

  • method ({‘nearest’, ‘linear’, ‘cubic’}) – The interpolation method.

  • extrapolation ({‘constant’, ‘nan’}) – The extrapolation method used for values of newx outside the range of xq. See notes.

Notes

Extrapolation methods:

  • ‘nan’ : Any value of newx outside the range of xq is set to NaN.

  • ‘constant’ : Values of newx smaller than the minimum of xq are set to the first value of yq and those larger than the maximum, set to the last one (first and last values along the “quantiles” dimension). When the grouping is “time.month”, these limits are linearly interpolated along the month dimension.

xclim.sdba.utils.invert(x: xarray.DataArray, kind: Optional[str] = None) xarray.DataArray[source]

Invert a DataArray either additively (-x) or multiplicatively (1/x).

If kind is not given, default to the one stored in the “kind” attribute of x.

xclim.sdba.utils.map_cdf(ds: xarray.Dataset, *, y_value: xarray.DataArray, dim)[source]

Return the value in x with the same CDF as y_value in y.

This function is meant to be wrapped in a Grouper.apply.

Parameters
  • ds (xr.Dataset) – Variables: x, Values from which to pick, y, Reference values giving the ranking

  • y_value (float, array) – Value within the support of y.

  • dim (str) – Dimension along which to compute quantile.

Returns

array – Quantile of x with the same CDF as y_value in y.

xclim.sdba.utils.map_cdf_1d(x, y, y_value)[source]

Return the value in x with the same CDF as y_value in y.

xclim.sdba.utils.pc_matrix(arr: Union[numpy.ndarray, dask.array.core.Array]) Union[numpy.ndarray, dask.array.core.Array][source]

Construct a Principal Component matrix.

This matrix can be used to transform points in arr to principal components coordinates. Note that this function does not manage NaNs; if a single observation is null, all elements of the transformation matrix involving that variable will be NaN.

Parameters

arr (numpy.ndarray or dask.array.Array) – 2D array (M, N) of the M coordinates of N points.

Returns

numpy.ndarray or dask.array.Array – MxM Array of the same type as arr.

xclim.sdba.utils.rand_rot_matrix(crd: xarray.DataArray, num: int = 1, new_dim: Optional[str] = None) xarray.DataArray[source]

Generate random rotation matrices.

Rotation matrices are members of the SO(n) group, where n is the matrix size (crd.size). They can be characterized as orthogonal matrices with determinant 1. A square matrix \(R\) is a rotation matrix if and only if \(R^t = R^{−1}\) and \(\mathrm{det} R = 1\).

Parameters
  • crd (xr.DataArray) – 1D coordinate DataArray along which the rotation occurs. The output will be square with the same coordinate replicated, the second renamed to new_dim.

  • num (int) – If larger than 1 (default), the number of matrices to generate, stacked along a “matrices” dimension.

  • new_dim (str) – Name of the new “prime” dimension, defaults to the same name as crd + “_prime”.

Returns

xr.DataArray – float, NxN if num = 1, numxNxN otherwise, where N is the length of crd.

References

Mezzadri, F. (2006). How to generate random matrices from the classical compact groups. arXiv preprint math-ph/0609050.

xclim.sdba.utils.rank(da: xarray.DataArray, dim: str = 'time', pct: bool = False) xarray.DataArray[source]

Ranks data along a dimension.

Replicates xr.DataArray.rank but as a function usable in a Grouper.apply(). Xarray’s docstring is below:

Equal values are assigned a rank that is the average of the ranks that would have been otherwise assigned to all the values within that set. Ranks begin at 1, not 0. If pct, computes percentage ranks.

Parameters
  • da (xr.DataArray) – Source array.

  • dim (str, hashable) – Dimension over which to compute rank.

  • pct (bool, optional) – If True, compute percentage ranks, otherwise compute integer ranks.

Returns

DataArray – DataArray with the same coordinates and dtype ‘float64’.

Notes

The bottleneck library is required. NaNs in the input array are returned as NaNs.

class xclim.sdba.base.Grouper(group: str, window: int = 1, add_dims: Optional[Union[Sequence[str], Set[str]]] = None)[source]

Create the Grouper object.

Parameters
  • group (str) – The usual grouping name as xarray understands it. Ex: “time.month” or “time”. The dimension name before the dot is the “main dimension” stored in Grouper.dim and the property name after is stored in Grouper.prop.

  • window (int) – If larger than 1, a centered rolling window along the main dimension is created when grouping data. Units are the sampling frequency of the data along the main dimension.

  • add_dims (Optional[Union[Sequence[str], str]]) – Additional dimensions that should be reduced in grouping operations. This behaviour is also controlled by the main_only parameter of the apply method. If any of these dimensions are absent from the dataarrays, they will be omitted.

apply(func: Union[function, str], da: Union[xarray.DataArray, Mapping[str, xarray.DataArray], xarray.Dataset], main_only: bool = False, **kwargs)[source]

Apply a function group-wise on DataArrays.

Parameters
  • func (Union[FunctionType, str]) – The function to apply to the groups, either a callable or a xr.core.groupby.GroupBy method name as a string. The function will be called as func(group, dim=dims, **kwargs). See main_only for the behaviour of dims.

  • da (Union[xr.DataArray, Mapping[str, xr.DataArray], xr.Dataset]) – The DataArray on which to apply the function. Multiple arrays can be passed through a dictionary. A dataset will be created before grouping.

  • main_only (bool) – Whether to call the function with the main dimension only (if True) or with all grouping dims (if False, default) (including the window and dimensions given through add_dims). The dimensions used are also written in the “group_compute_dims” attribute. If all the input arrays are missing one of the ‘add_dims’, it is silently omitted.

  • kwargs – Other keyword arguments to pass to the function.

Returns

DataArray or Dataset – Attributes “group”, “group_window” and “group_compute_dims” are added.

If the function did not reduce the array:

  • The output is sorted along the main dimension.

  • The output is rechunked to match the chunks on the input If multiple inputs with differing chunking were given as inputs, the chunking with the smallest number of chunks is used.

If the function reduces the array:

  • If there is only one group, the singleton dimension is squeezed out of the output

  • The output is rechunked as to have only 1 chunk along the new dimension.

Notes

For the special case where a Dataset is returned, but only some of its variable where reduced by the grouping, xarray’s GroupBy.map will broadcast everything back to the ungrouped dimensions. To overcome this issue, function may add a “_group_apply_reshape” attribute set to True on the variables that should be reduced and these will be re-grouped by calling da.groupby(self.name).first().

property freq

The frequency string corresponding to the group. For use with xarray’s resampling.

get_coordinate(ds=None)[source]

Return the coordinate as in the output of group.apply.

Currently, only implemented for groupings with prop == month or dayofyear. For prop == dayfofyear, a ds (Dataset or DataArray) can be passed to infer the max doy from the available years and calendar.

get_index(da: Union[xarray.DataArray, xarray.Dataset], interp: Optional[bool] = None)[source]

Return the group index of each element along the main dimension.

Parameters
  • da (Union[xr.DataArray, xr.Dataset]) – The input array/dataset for which the group index is returned. It must have Grouper.dim as a coordinate.

  • interp (bool, optional) – If True, the returned index can be used for interpolation. Only value for month grouping, where integer values represent the middle of the month, all other days are linearly interpolated in between.

Returns

xr.DataArray – The index of each element along Grouper.dim. If Grouper.dim is time and Grouper.prop is None, an uniform array of True is returned. If Grouper.prop is a time accessor (month, dayofyear, etc), an numerical array is returned, with a special case of month and interp=True. If Grouper.dim is not time, the dim is simply returned.

group(da: Optional[Union[xarray.DataArray, xarray.Dataset]] = None, main_only=False, **das: xarray.DataArray)[source]

Return a xr.core.groupby.GroupBy object.

More than one array can be combined to a dataset before grouping using the das kwargs. A new window dimension is added if self.window is larger than 1. If Grouper.dim is ‘time’, but ‘prop’ is None, the whole array is grouped together.

When multiple arrays are passed, some of them can be grouped along the same group as self. They are boadcasted, merged to the grouping dataset and regrouped in the output.

property prop_name

A significative name of the grouping.

Numba-accelerated utilities

xclim.sdba.nbutils.quantile(da, q, dim)[source]

Compute the quantiles from a fixed list “q”

xclim.sdba.nbutils.vecquantiles(da, rnk, dim)[source]

For when the quantile (rnk) is different for each point.

da and rnk must share all dimensions but dim.

LOESS smoothing

xclim.sdba.loess.loess_smoothing(da: xarray.DataArray, dim: str = 'time', d: int = 1, f: float = 0.5, niter: int = 2, weights: Union[str, Callable] = 'tricube', equal_spacing: Optional[bool] = None)[source]

Locally weighted regression in 1D: fits a nonparametric regression curve to a scatterplot.

Returns a smoothed curve along given dimension. The regression is computed for each point using a subset of neigboring points as given from evaluating the weighting function locally. Follows the procedure of [Cleveland1979].

Parameters
  • da (xr.DataArray) – The data to smooth using the loess approach.

  • dim (str) – Name of the dimension along which to perform the loess.

  • d ([0, 1]) – Degree of the local regression.

  • f (float) – Parameter controlling the shape of the weight curve. Behavior depends on the weighting function, but it usually represents the span of the weighting function in reference to x-coordinates normalized from 0 to 1.

  • niter (int) – Number of robustness iterations to execute.

  • weights ([“tricube”, “gaussian”] or callable) – Shape of the weighting function, see notes. The user can provide a function or a string: “tricube” : a smooth top-hat like curve. “gaussian” : a gaussian curve, f gives the span for 95% of the values.

  • equal_spacing (bool, optional) – Whether to use the equal spacing optimization. If None (the default), it is activated only if the x-axis is equally-spaced. When activated, dx = x[1] - x[0].

Notes

As stated in [Cleveland1979], the weighting function \(W(x)\) should respect the following conditions:

  • \(W(x) > 0\) for \(|x| < 1\)

  • \(W(-x) = W(x)\)

  • \(W(x)\) is non-increasing for \(x \ge 0\)

  • \(W(x) = 0\) for \(|x| \ge 0\)

If a Callable is provided, it should only accept the 1D np.ndarray \(x\) which is an absolute value function going from 1 to 0 to 1 around \(x_i\), for all values where \(x - x_i < h_i\) with \(h_i\) the distance of the rth nearest neighbor of \(x_i\), \(r = f * size(x)\).

References

Cleveland1979(1,2)

Cleveland, W. S., 1979. Robust Locally Weighted Regression and Smoothing Scatterplot, Journal of the American Statistical Association 74, 829–836.

Code adapted from https://gist.github.com/agramfort/850437

Properties submodule

SDBA diagnostic tests are made up of statistical properties and measures. Properties are calculated on both simulation and reference datasets. They collapse the time dimension to one value.

This framework for the diagnostic tests was inspired by the [VALUE] project. Statistical Properties is the xclim term for ‘indices’ in the VALUE project.

VALUE

http://www.value-cost.eu/

xclim.sdba.properties.STATISTICAL_PROPERTIES: Dict[str, Callable] = {'acf': <function acf>, 'annual_cycle_amplitude': <function annual_cycle_amplitude>, 'annual_cycle_phase': <function annual_cycle_phase>, 'corr_btw_var': <function corr_btw_var>, 'mean': <function mean>, 'quantile': <function quantile>, 'relative_frequency': <function relative_frequency>, 'return_value': <function return_value>, 'skewness': <function skewness>, 'spell_length_distribution': <function spell_length_distribution>, 'trend': <function trend>, 'var': <function var>}

Dictionary of all the statistical properties available.

xclim.sdba.properties.acf(da: xarray.DataArray, *, lag: int = 1, group: Union[str, xclim.sdba.base.Grouper] = 'time.season') xarray.DataArray[source]

Autocorrelation function.

Autocorrelation with a lag over a time resolution and averaged over all years.

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • lag (int) – Lag.

  • group ({‘time.season’, ‘time.month’}) – Grouping of the output. E.g. If ‘time.month’, the autocorrelation is calculated over each month separately for all years. Then, the autocorrelation for all Jan/Feb/… is averaged over all years, giving 12 outputs for each grid point.

Returns

xr.DataArray – lag-{lag} autocorrelation of the variable over a {group.prop} and averaged over all years.

References

Alavoine M., and Grenier P. (under review) The distinct problems of physical inconsistency and of multivariate bias potentially involved in the statistical adjustment of climate simulations. International Journal of Climatology, submitted on September 19th 2021. (Preprint: https://doi.org/10.31223/X5C34C)

Examples

>>> from xclim.testing import open_dataset
>>> pr = open_dataset(path_to_pr_file).pr
>>> acf(da=pr, lag=3, group='time.season')
xclim.sdba.properties.annual_cycle_amplitude(da: xarray.DataArray, *, amplitude_type: str = 'absolute', group: Union[str, xclim.sdba.base.Grouper] = 'time') xarray.DataArray[source]

Annual cycle amplitude.

The amplitudes of the annual cycle are calculated for each year, than averaged over the all years.

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • amplitude_type ({‘absolute’,’relative’}) – Type of amplitude. ‘absolute’ is the peak-to-peak amplitude. (max - min). ‘relative’ is a relative percentage. 100 * (max - min) / mean (Recommended for precipitation).

Returns

xr.DataArray – {amplitude_type} amplitude of the annual cycle.

Examples

>>> pr = open_dataset(path_to_pr_file).pr
>>> annual_cycle_amplitude(da=pr, amplitude_type='relative')
xclim.sdba.properties.annual_cycle_phase(da: xarray.DataArray, *, group: Union[str, xclim.sdba.base.Grouper] = 'time') xarray.DataArray[source]

Annual cycle phase.

The phases of the annual cycle are calculated for each year, than averaged over the all years.

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • group ({“time”, ‘time.season’, ‘time.month’}) – Grouping of the output. Default: “time”.

Returns

xr.DataArray – Phase of the annual cycle. The position (day-of-year) of the maximal value.

Examples

>>> pr = open_dataset(path_to_pr_file).pr
>>> annual_cycle_phase(da=pr)
xclim.sdba.properties.corr_btw_var(da1: xarray.DataArray, da2: xarray.DataArray, *, corr_type: str = 'Spearman', group: Union[str, xclim.sdba.base.Grouper] = 'time', output: str = 'correlation') xarray.DataArray[source]

Correlation between two variables.

Spearman or Pearson correlation coefficient between two variables at the time resolution.

Parameters
  • da1 (xr.DataArray) – First variable on which to calculate the diagnostic.

  • da2 (xr.DataArray) – Second variable on which to calculate the diagnostic.

  • corr_type ({‘Pearson’,’Spearman’}) – Type of correlation to calculate.

  • output ({‘correlation’, ‘pvalue’}) – Wheter to return the correlation coefficient or the p-value.

  • group ({‘time’, ‘time.season’, ‘time.month’}) – Grouping of the output. Eg. For ‘time.month’, the correlation would be calculated on each month separately, but with all the years together.

Returns

xr.DataArray – {corr_type} correlation coefficient

Examples

>>> pr = open_dataset(path_to_pr_file).pr
>>> tasmax = open_dataset('NRCANdaily/nrcan_canada_daily_tasmax_1990.nc').tasmax
>>> corr_btw_var(da1=pr, da2=tasmax, group='time.season')
xclim.sdba.properties.mean(da: xarray.DataArray, *, group: Union[str, xclim.sdba.base.Grouper] = 'time') xarray.DataArray[source]

Mean.

Mean over all years at the time resolution.

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • group ({‘time’, ‘time.season’, ‘time.month’}) – Grouping of the output. E.g. If ‘time.month’, the temporal average is performed separately for each month.

Returns

xr.DataArray, – Mean of the variable.

Examples

>>> pr = open_dataset(path_to_pr_file).pr
>>> mean(da=pr, group='time.season')
xclim.sdba.properties.quantile(da: xarray.DataArray, *, q: float = 0.98, group: Union[str, xclim.sdba.base.Grouper] = 'time') xarray.DataArray[source]

Quantile.

Returns the quantile q of the distribution of the variable over all years at the time resolution.

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • q (float) – Quantile to be calculated. Should be between 0 and 1.

  • group ({‘time’, ‘time.season’, ‘time.month’}) – Grouping of the output. E.g. If ‘time.month’, the quantile is computed separately for each month.

Returns

xr.DataArray – Quantile {q} of the variable.

Examples

>>> pr = open_dataset(path_to_pr_file).pr
>>> quantile(da=pr, q=0.9, group='time.season')
xclim.sdba.properties.relative_frequency(da: xarray.DataArray, *, op: str = '>=', thresh: str = '1mm d-1', group: Union[str, xclim.sdba.base.Grouper] = 'time') xarray.DataArray[source]

Relative Frequency.

Relative Frequency of days with variable respecting a condition (defined by an operation and a threshold) at the time resolution. The relative freqency is the number of days that satisfy the condition divided by the total number of days.

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • op ({“>”, “<”, “>=”, “<=”}) – Operation to verify the condition. The condition is variable {op} threshold.

  • thresh (str) – Threshold on which to evaluate the condition.

  • group ({‘time’, ‘time.season’, ‘time.month’}) – Grouping on the output. Eg. For ‘time.month’, the relative frequency would be calculated on each month, with all years included.

Returns

xr.DataArray – Relative frequency of the variable.

Examples

>>> tasmax = open_dataset(path_to_tasmax_file).tasmax
>>> relative_frequency(da=tasmax, op= '<', thresh= '0 degC', group='time.season')
xclim.sdba.properties.return_value(da: xarray.DataArray, *, period: int = 20, op: str = 'max', method: str = 'ML', group: Union[str, xclim.sdba.base.Grouper] = 'time') xarray.DataArray[source]

Return value.

Return the value corresponding to a return period. On average, the return value will be exceeded (or not exceed for op=’min’) every return period (eg. 20 years). The return value is computed by first extracting the variable annual maxima/minima, fitting a statistical distribution to the maxima/minima, then estimating the percentile associated with the return period (eg. 95th percentile (1/20) for 20 years)

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • period (int) – Return period. Number of years over which to check if the value is exceeded (or not for op=’min’).

  • op ({‘max’,’min’}) – Whether we are looking for a probability of exceedance (‘max’, right side of the distribution) or a probability of non-exceedance (min, left side of the distribution).

  • method ({“ML”, “PWM”}) – Fitting method, either maximum likelihood (ML) or probability weighted moments (PWM), also called L-Moments. The PWM method is usually more robust to outliers. However, it requires the lmoments3 libraryto be installed from the develop branch. pip install git+https://github.com/OpenHydrology/lmoments3.git@develop#egg=lmoments3

  • group ({‘time’, ‘time.season’, ‘time.month’}) – Grouping of the output. A distribution of the extremums is done for each group.

Returns

xr.DataArray – {period}-{group} {op} return level of the variable.

Examples

>>> tas = open_dataset(path_to_tas_file).tas
>>> return_value(da=tas, group='time.season')
xclim.sdba.properties.skewness(da: xarray.DataArray, *, group: Union[str, xclim.sdba.base.Grouper] = 'time') xarray.DataArray[source]

Skewness.

Skewness of the distribution of the variable over all years at the time resolution.

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • group ({‘time’, ‘time.season’, ‘time.month’}) – Grouping of the output. E.g. If ‘time.month’, the skewness is performed separately for each month.

Returns

xr.DataArray – Skewness of the variable.

Examples

>>> pr = open_dataset(path_to_pr_file).pr
>>> skewness(da=pr, group='time.season')

See also

scipy.stats.skew

xclim.sdba.properties.spell_length_distribution(da: xarray.DataArray, *, method: str = 'amount', op: str = '>=', thresh: Union[str, float] = '1 mm d-1', stat: str = 'mean', group: Union[str, xclim.sdba.base.Grouper] = 'time') xarray.DataArray[source]

Spell length distribution.

Statistic of spell length distribution when the variable respects a condition (defined by an operation, a method and

a threshold).

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • method ({‘amount’, ‘quantile’}) – Method to choose the threshold. ‘amount’: The threshold is directly the quantity in {thresh}. It needs to have the same units as {da}. ‘quantile’: The threshold is calculated as the quantile {thresh} of the distribution.

  • op ({“>”, “<”, “>=”, “<=”}) – Operation to verify the condition for a spell. The condition for a spell is variable {op} threshold.

  • thresh (str or float) – Threshold on which to evaluate the condition to have a spell. Str with units if the method is “amount”. Float of the quantile if the method is “quantile”.

  • stat ({‘mean’,’max’,’min’}) – Statistics to apply to the resampled input at the {group} (e.g. 1-31 Jan 1980) and then over all years (e.g. Jan 1980-2010)

  • group ({‘time’, ‘time.season’, ‘time.month’}) – Grouping of the output. E.g. If ‘time.month’, the spell lengths are coputed separately for each month.

Returns

xr.DataArray – {stat} of spell length distribution when the variable is {op} the {method} {thresh}.

Examples

>>> pr = open_dataset(path_to_pr_file).pr
>>> spell_length_distribution(da=pr, op='<',thresh ='1mm d-1', group='time.season')
xclim.sdba.properties.trend(da: xarray.DataArray, *, group: Union[str, xclim.sdba.base.Grouper] = 'time', output: str = 'slope') xarray.DataArray[source]

Linear Trend.

The data is averaged over each time resolution and the interannual trend is returned.

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • output ({‘slope’, ‘pvalue’}) – Attributes of the linear regression to return. ‘slope’ is the slope of the regression line. ‘pvalue’ is for a hypothesis test whose null hypothesis is that the slope is zero, using Wald Test with t-distribution of the test statistic.

  • group ({‘time’, ‘time.season’, ‘time.month’}) – Grouping on the output.

Returns

xr.DataArray – Trend of the variable.

Examples

>>> tas = open_dataset(path_to_tas_file).tas
>>> trend(da=tas, group='time.season')
xclim.sdba.properties.var(da: xarray.DataArray, *, group: Union[str, xclim.sdba.base.Grouper] = 'time') xarray.DataArray[source]

Variance.

Variance of the variable over all years at the time resolution.

Parameters
  • da (xr.DataArray) – Variable on which to calculate the diagnostic.

  • group ({‘time’, ‘time.season’, ‘time.month’}) – Grouping of the output. E.g. If ‘time.month’, the variance is performed separately for each month.

Returns

xr.DataArray – Variance of the variable.

Examples

>>> pr = open_dataset(path_to_pr_file).pr
>>> var(da=pr, group='time.season')

Measures submodule

SDBA diagnostic tests are made up of properties and measures. Measures compare adjusted simulations to a reference, through statistical properties or directly. This framework for the diagnostic tests was inspired by the [VALUE] project.

xclim.sdba.measures.annual_cycle_correlation(sim, ref, window: int = 15)[source]

Annual cycle correlation.

Pearson correlation coefficient between the smooth day-of-year averaged annual cycles of the simulation and the reference. In the smooth day-of-year averaged annual cycles, each day-of-year is averaged over all years and over a window of days around that day.

Parameters
  • sim (xr.DataArray) – data from the simulation (a time-series for each grid-point)

  • ref (xr.DataArray) – data from the reference (observations) (a time-series for each grid-point)

  • window (int) – Size of window around each day of year around which to take the mean. E.g. If window=31, Jan 1st is averaged over from December 17th to January 16th.

Returns

xr.DataArray, – Annual cycle correlation between the simulation and the reference

xclim.sdba.measures.bias(sim: xarray.DataArray, ref: xarray.DataArray) xarray.DataArray[source]

Bias.

The bias is the simulation minus the reference.

Parameters
  • sim (xr.DataArray) – data from the simulation (one value for each grid-point)

  • ref (xr.DataArray) – data from the reference (observations) (one value for each grid-point)

Returns

xr.DataArray, – Bias between the simulation and the reference

xclim.sdba.measures.circular_bias(sim: xarray.DataArray, ref: xarray.DataArray) xarray.DataArray[source]

Circular bias.

Bias considering circular time series. E.g. The bias between doy 365 and doy 1 is 364, but the circular bias is -1.

Parameters
  • sim (xr.DataArray) – data from the simulation (one value for each grid-point)

  • ref (xr.DataArray) – data from the reference (observations) (one value for each grid-point)

Returns

xr.DataArray, – Circular bias between the simulation and the reference

xclim.sdba.measures.mae(sim: xarray.DataArray, ref: xarray.DataArray) xarray.DataArray[source]

Mean absolute error.

The mean absolute error on the time dimension between the simulation and the reference.

Parameters
  • sim (xr.DataArray) – data from the simulation (a time-series for each grid-point)

  • ref (xr.DataArray) – data from the reference (observations) (a time-series for each grid-point)

Returns

xr.DataArray, – Mean absolute error between the simulation and the reference

xclim.sdba.measures.ratio(sim: xarray.DataArray, ref: xarray.DataArray) xarray.DataArray[source]

Ratio.

The ratio is the quotient of the simulation over the reference.

Parameters
  • sim (xr.DataArray) – data from the simulation (one value for each grid-point)

  • ref (xr.DataArray) – data from the reference (observations) (one value for each grid-point)

Returns

xr.DataArray, – Ratio between the simulation and the reference

xclim.sdba.measures.relative_bias(sim: xarray.DataArray, ref: xarray.DataArray) xarray.DataArray[source]

Relative Bias.

The relative bias is the simulation minus reference, divided by the reference.

Parameters
  • sim (xr.DataArray) – data from the simulation (one value for each grid-point)

  • ref (xr.DataArray) – data from the reference (observations) (one value for each grid-point)

Returns

xr.DataArray, – Relative bias between the simulation and the reference

xclim.sdba.measures.rmse(sim: xarray.DataArray, ref: xarray.DataArray) xarray.DataArray[source]

Root mean square error.

The root mean square error on the time dimension between the simulation and the reference.

Parameters
  • sim (xr.DataArray) – Data from the simulation (a time-series for each grid-point)

  • ref (xr.DataArray) – Data from the reference (observations) (a time-series for each grid-point)

Returns

xr.DataArray, – Root mean square error between the simulation and the reference

Developer tools

Base classes and developper tools

class xclim.sdba.base.Parametrizable[source]

Bases: dict

Helper base class resembling a dictionary.

This object is _completely_ defined by the content of its internal dictionary, accessible through item access (self[‘attr’]) or in self.parameters. When serializing and restoring this object, only members of that internal dict are preserved. All other attributes set directly with self.attr = value will not be preserved upon serialization and restoration of the object with [json]pickle. dictionary. Other variables set with self.var = data will be lost in the serialization process. This class is best serialized and restored with jsonpickle.

property parameters

All parameters as a dictionary. Read-only.

class xclim.sdba.base.ParametrizableWithDataset[source]

Bases: xclim.sdba.base.Parametrizable

Parametrizeable class that also has a ds attribute storing a dataset.

classmethod from_dataset(ds: xarray.Dataset)[source]

Create an instance from a dataset.

The dataset must have a global attribute with a name corresponding to cls._attribute, and that attribute must be the result of jsonpickle.encode(object) where object is of the same type as this object.

set_dataset(ds: xarray.Dataset)[source]

Stores an xarray dataset in the ds attribute.

Useful with custom object initialization or if some external processing was performed.

xclim.sdba.base.duck_empty(dims, sizes, dtype='float64', chunks=None)[source]

Return an empty DataArray based on a numpy or dask backend, depending on the chunks argument.

xclim.sdba.base.map_blocks(reduces: Optional[Sequence[str]] = None, **outvars)[source]

Decorator for declaring functions and wrapping them into a map_blocks. It takes care of constructing the template dataset.

Dimension order is not preserved.

The decorated function must always have the signature: func(ds, **kwargs), where ds is a DataArray or a Dataset. It must always output a dataset matching the mapping passed to the decorator.

Parameters
  • reduces (sequence of strings) – Name of the dimensions that are removed by the function.

  • outvars – Mapping from variable names in the output to their new dimensions. The placeholders Grouper.PROP, Grouper.DIM and Grouper.ADD_DIMS can be used to signify group.prop,`group.dim` and group.add_dims respectively. If an output keeps a dimension that another loses, that dimension name must be given in reduces and in the list of new dimensions of the first output.

xclim.sdba.base.map_groups(reduces: Optional[Sequence[str]] = None, main_only: bool = False, **outvars)[source]

Decorator for declaring functions acting only on groups and wrapping them into a map_blocks. See map_blocks().

This is the same as map_blocks but adds a call to group.apply() in the mapped func and the default value of reduces is changed.

The decorated function must have the signature: func(ds, dim, **kwargs). Where ds is a DataAray or Dataset, dim is the group.dim (and add_dims). The group argument is stripped from the kwargs, but must evidently be provided in the call.

Parameters
  • reduces (sequence of str) – Dimensions that are removed from the inputs by the function. Defaults to [Grouper.DIM, Grouper.ADD_DIMS] if main_only is False, and [Grouper.DIM] if main_only is True. See map_blocks().

  • main_only (bool) – Same as for Grouper.apply().

  • outvars – Mapping from variable names in the output to their new dimensions. The placeholders Grouper.PROP, Grouper.DIM and Grouper.ADD_DIMS can be used to signify group.prop,`group.dim` and group.add_dims respectively. If an output keeps a dimension that another loses, that dimension name must be given in reduces and in the list of new dimensions of the first output.

xclim.sdba.base.parse_group(func: Callable, kwargs=None, allow_only=None) Callable[source]

Parse the kwargs given to a function to set the group arg with a Grouper object.

This function can be used as a decorator, in which case the parsing and updating of the kwargs is done at call time. It can also be called with a function from which extract the default group and kwargs to update, in which case it returns the updated kwargs.

If allow_only is given, an exception is raised when the parsed group is not within that list.

class xclim.sdba.detrending.BaseDetrend(*, group: Union[xclim.sdba.base.Grouper, str] = 'time', kind: str = '+', **kwargs)[source]

Base class for detrending objects.

Defines three methods:

fit(da) : Compute trend from da and return a new _fitted_ Detrend object. detrend(da) : Return detrended array. retrend(da) : Puts trend back on da.

A fitted Detrend object is unique to the trend coordinate of the object used in fit, (usually ‘time’). The computed trend is stored in Detrend.ds.trend.

Subclasses should implement _get_trend_group() or _get_trend(). The first will be called in a group.apply(..., main_only=True), and should return a single DataArray. The second allows the use of functions wrapped in map_groups() and should also return a single DataArray.

The subclasses may reimplement _detrend and _retrend.

detrend(da: xarray.DataArray)[source]

Remove the previously fitted trend from a DataArray.

fit(da: xarray.DataArray)[source]

Extract the trend of a DataArray along a specific dimension.

Returns a new object that can be used for detrending and retrending. Fitted objects are unique to the fitted coordinate used.

retrend(da: xarray.DataArray)[source]

Put the previously fitted trend back on a DataArray.

class xclim.sdba.adjustment.TrainAdjust(*args, _trained=False, **kwargs)[source]

Base class for adjustment objects obeying the train-adjust scheme.

Children classes should implement these methods:

  • _train(ref, hist, **kwargs), classmethod receiving the training target and data, returning a training dataset and parameters to store in the object.

  • _adjust(sim, **kwargs), receiving the projected data and some arguments, returning the scen dataarray.

adjust(sim: xarray.DataArray, *args, **kwargs)[source]

Return bias-adjusted data. Refer to the class documentation for the algorithm details.

Parameters
  • sim (DataArray) – Time series to be bias-adjusted, usually a model output.

  • args (xr.DataArray) – Other DataArrays needed for the adjustment (usually none).

  • kwargs – Algorithm-specific keyword arguments, see class doc.

set_dataset(ds: xarray.Dataset)[source]

Stores an xarray dataset in the ds attribute.

Useful with custom object initialization or if some external processing was performed.

classmethod train(ref: xarray.DataArray, hist: xarray.DataArray, **kwargs)[source]

Train the adjustment object. Refer to the class documentation for the algorithm details.

Parameters
  • ref (DataArray) – Training target, usually a reference time series drawn from observations.

  • hist (DataArray) – Training data, usually a model output whose biases are to be adjusted.

class xclim.sdba.adjustment.Adjust(*args, _trained=False, **kwargs)[source]

Adjustment with no intermediate trained object.

Children classes should implement a _adjust classmethod taking as input the three DataArrays and returning the scen dataset/array.

classmethod adjust(ref: xarray.DataArray, hist: xarray.DataArray, sim: xarray.DataArray, **kwargs)[source]

Return bias-adjusted data. Refer to the class documentation for the algorithm details.

Parameters
  • ref (DataArray) – Training target, usually a reference time series drawn from observations.

  • hist (DataArray) – Training data, usually a model output whose biases are to be adjusted.

  • sim (DataArray) – Time series to be bias-adjusted, usually a model output.

  • kwargs – Algorithm-specific keyword arguments, see class doc.

xclim.sdba.properties.register_statistical_properties(aspect: str, seasonal: bool, annual: bool) Callable[source]

Register statistical properties in the STATISTICAL_PROPERTIES dictionary with its aspect and time resolutions.

xclim.sdba.measures.check_same_units_and_convert(func) Callable[source]

Verify that the simulation and the reference have the same units. If not, it converts the simulation to the units of the reference