xclim package

Climate indices computation package based on Xarray.

Subpackages

Submodules

xclim.analog module

Spatial Analogues module.

xclim.analog.friedman_rafsky(x, y)[source]

Compute a dissimilarity metric based on the Friedman-Rafsky runs statistics.

The algorithm builds a minimal spanning tree (the subset of edges connecting all points that minimizes the total edge length) then counts the edges linking points from the same distribution. This method is scale-dependent.

Parameters:
  • x (np.ndarray (n,d)) – Reference sample.

  • y (np.ndarray (m,d)) – Candidate sample.

Return type:

float

Returns:

float – Friedman-Rafsky dissimilarity metric ranging from 0 to (m+n-1)/(m+n).

References

Friedman and Rafsky [1979]

xclim.analog.kldiv(x, y, *, k=1)[source]

Compute the Kullback-Leibler divergence between two multivariate samples.

The formula to compute the K-L divergence from samples is given by:

\[D(P||Q) = \frac{d}{n} \sum_i^n \log\left\{\frac{r_k(x_i)}{s_k(x_i)}\right\} + \log\left\{\frac{m}{n-1}\right\}\]

where \(r_k(x_i)\) and \(s_k(x_i)\) are, respectively, the Euclidean distance to the kth neighbour of \(x_i\) in the x array (excepting \(x_i\)) and in the y array. This method is scale-dependent.

Parameters:
  • x (np.ndarray (n,d)) – Samples from distribution P, which typically represents the true distribution (reference).

  • y (np.ndarray (m,d)) – Samples from distribution Q, which typically represents the approximate distribution (candidate)

  • k (int or sequence) – The kth neighbours to look for when estimating the density of the distributions. Defaults to 1, which can be noisy.

Return type:

float | Sequence[float]

Returns:

float or sequence – The estimated Kullback-Leibler divergence D(P||Q) computed from the distances to the kth neighbour.

Notes

In information theory, the Kullback–Leibler divergence [Perez-Cruz, 2008] is a non-symmetric measure of the difference between two probability distributions P and Q, where P is the “true” distribution and Q an approximation. This nuance is important because \(D(P||Q)\) is not equal to \(D(Q||P)\).

For probability distributions P and Q of a continuous random variable, the K–L divergence is defined as:

\[D_{KL}(P||Q) = \int p(x) \log\left(\frac{p(x)}{q(x)}\right) dx\]

This formula assumes we have a representation of the probability densities \(p(x)\) and \(q(x)\). In many cases, we only have samples from the distribution, and most methods first estimate the densities from the samples and then proceed to compute the K-L divergence. In Perez-Cruz [2008], the author proposes an algorithm to estimate the K-L divergence directly from the sample using an empirical CDF. Even though the CDFs do not converge to their true values, the paper proves that the K-L divergence almost surely does converge to its true value.

References

Perez-Cruz [2008]

xclim.analog.kolmogorov_smirnov(x, y)[source]

Compute the Kolmogorov-Smirnov statistic applied to two multivariate samples as described by Fasano and Franceschini.

This method is scale-dependent.

Parameters:
  • x (np.ndarray (n,d)) – Reference sample.

  • y (np.ndarray (m,d)) – Candidate sample.

Return type:

float

Returns:

float – Kolmogorov-Smirnov dissimilarity metric ranging from 0 to 1.

References

Fasano and Franceschini [1987]

xclim.analog.metric(func)[source]

Register a metric function in the metrics mapping and add some preparation/checking code.

All metric functions accept 2D inputs. This reshapes 1D inputs to (n, 1) and (m, 1). All metric functions are invalid when any non-finite values are present in the inputs.

xclim.analog.nearest_neighbor(x, y)[source]

Compute a dissimilarity metric based on the number of points in the pooled sample whose nearest neighbor belongs to the same distribution.

This method is scale-invariant.

Parameters:
  • x (np.ndarray (n,d)) – Reference sample.

  • y (np.ndarray (m,d)) – Candidate sample.

Return type:

ndarray

Returns:

float – Nearest-Neighbor dissimilarity metric ranging from 0 to 1.

References

Henze [1988]

xclim.analog.seuclidean(x, y)[source]

Compute the Euclidean distance between the mean of a multivariate candidate sample with respect to the mean of a reference sample.

This method is scale-invariant.

Parameters:
  • x (np.ndarray (n,d)) – Reference sample.

  • y (np.ndarray (m,d)) – Candidate sample.

Return type:

float

Returns:

float – Standardized Euclidean Distance between the mean of the samples ranging from 0 to infinity.

Notes

This metric considers neither the information from individual points nor the standard deviation of the candidate distribution.

References

Veloz, Williams, Lorenz, Notaro, Vavrus, and Vimont [2012]

xclim.analog.spatial_analogs(target, candidates, dist_dim='time', method='kldiv', **kwargs)[source]

Compute dissimilarity statistics between target points and candidate points.

Spatial analogues based on the comparison of climate indices. The algorithm compares the distribution of the reference indices with the distribution of spatially distributed candidate indices and returns a value measuring the dissimilarity between both distributions over the candidate grid.

Parameters:
  • target (xr.Dataset) – Dataset of the target indices. Only indice variables should be included in the dataset’s data_vars. They should have only the dimension(s) dist_dim `in common with `candidates.

  • candidates (xr.Dataset) – Dataset of the candidate indices. Only indice variables should be included in the dataset’s data_vars.

  • dist_dim (str) – The dimension over which the distributions are constructed. This can be a multi-index dimension.

  • method ({‘seuclidean’, ‘nearest_neighbor’, ‘zech_aslan’, ‘kolmogorov_smirnov’, ‘friedman_rafsky’, ‘kldiv’}) – Which method to use when computing the dissimilarity statistic.

  • **kwargs – Any other parameter passed directly to the dissimilarity method.

Returns:

xr.DataArray – The dissimilarity statistic over the union of candidates’ and target’s dimensions. The range depends on the method.

xclim.analog.standardize(x, y)[source]

Standardize x and y by the square root of the product of their standard deviation.

Parameters:
  • x (np.ndarray) – Array to be compared.

  • y (np.ndarray) – Array to be compared.

Return type:

tuple[ndarray, ndarray]

Returns:

(ndarray, ndarray) – Standardized arrays.

xclim.analog.szekely_rizzo(x, y, *, standardize=True)[source]

Compute the Székely-Rizzo energy distance dissimilarity metric based on an analogy with Newton’s gravitational potential energy.

This method is scale-invariant when standardize=True (default), scale-dependent otherwise.

Parameters:
  • x (ndarray (n,d)) – Reference sample.

  • y (ndarray (m,d)) – Candidate sample.

  • standardize (bool) – If True (default), the standardized euclidean norm is used, instead of the conventional one.

Return type:

float

Returns:

float – Székely-Rizzo’s energy distance dissimilarity metric ranging from 0 to infinity.

Notes

The e-distance between two variables \(X\), \(Y\) (target and candidates) of sizes \(n,d\) and \(m,d\) proposed by Szekely and Rizzo [2004] is defined by:

\[e(X, Y) = \frac{n m}{n + m} \left[2\phi_{xy} − \phi_{xx} − \phi_{yy} \right]\]

where

\[\begin{split}\phi_{xy} &= \frac{1}{n m} \sum_{i = 1}^n \sum_{j = 1}^m \left\Vert X_i − Y_j \right\Vert \\ \phi_{xx} &= \frac{1}{n^2} \sum_{i = 1}^n \sum_{j = 1}^n \left\Vert X_i − X_j \right\Vert \\ \phi_{yy} &= \frac{1}{m^2} \sum_{i = 1}^m \sum_{j = 1}^m \left\Vert X_i − Y_j \right\Vert \\\end{split}\]

and where \(\Vert\cdot\Vert\) denotes the Euclidean norm, \(X_i\) denotes the i-th observation of \(X\). When standardized=False, this corresponds to the \(T\) test of Rizzo and Székely [2016] (p. 28) and to the eqdist.e function of the energy R package (with two samples) and gives results twice as big as xclim.sdba.processing.escore(). The standardization was added following the logic of [Grenier et al., 2013] to make the metric scale-invariant.

References

Grenier, Parent, Huard, Anctil, and Chaumont [2013], Rizzo and Székely [2016], Szekely and Rizzo [2004]

xclim.analog.zech_aslan(x, y, *, dmin=1e-12)[source]

Compute a modified Zech-Aslan energy distance dissimilarity metric based on an analogy with the energy of a cloud of electrical charges.

This method is scale-invariant.

Parameters:
  • x (np.ndarray (n,d)) – Reference sample.

  • y (np.ndarray (m,d)) – Candidate sample.

  • dmin (float) – The cut-off for low distances to avoid singularities on identical points.

Return type:

float

Returns:

float – Zech-Aslan dissimilarity metric ranging from -infinity to infinity.

Notes

The energy measure between two variables \(X\), \(Y\) (target and candidates) of sizes \(n,d\) and \(m,d\) proposed by Aslan and Zech [2003] is defined by:

\[\begin{split}e(X, Y) &= \left[\phi_{xx} + \phi_{yy} - \phi_{xy}\right] \\ \phi_{xy} &= \frac{1}{n m} \sum_{i = 1}^n \sum_{j = 1}^m R\left[SED(X_i, Y_j)\right] \\ \phi_{xx} &= \frac{1}{n^2} \sum_{i = 1}^n \sum_{j = i + 1}^n R\left[SED(X_i, X_j)\right] \\ \phi_{yy} &= \frac{1}{m^2} \sum_{i = 1}^m \sum_{j = i + 1}^m R\left[SED(X_i, Y_j)\right] \\\end{split}\]

where \(X_i\) denotes the i-th observation of \(X\). \(R\) is a weight function and \(SED(A, B)\) denotes the standardized Euclidean distance.

\[\begin{split}R(r) &= \left\{\begin{array}{r l} -\ln r & \text{for } r > d_{min} \\ -\ln d_{min} & \text{for } r \leq d_{min} \end{array}\right. \\ SED(X_i, Y_j) &= \sqrt{\sum_{k=1}^d \frac{\left(X_i(k) - Y_i(k)\right)^2}{\sigma_x(k)\sigma_y(k)}}\end{split}\]

where \(k\) is a counter over dimensions (indices in the case of spatial analogs) and \(\sigma_x(k)\) is the standard deviation of \(X\) in dimension \(k\). Finally, \(d_{min}\) is a cut-off to avoid poles when \(r \to 0\), it is controllable through the dmin parameter.

This version corresponds the \(D_{ZAE}\) test of Grenier et al. [2013] (eq. 7), which is a version of \(\phi_{NM}\) from Aslan and Zech [2003], modified by using the standardized Euclidean distance, the log weight function and choosing \(d_{min} = 10^{-12}\).

References

Aslan and Zech [2003], Grenier, Parent, Huard, Anctil, and Chaumont [2013], Zech and Aslan [2003]

xclim.cli module

Command Line Interface module

class xclim.cli.XclimCli(name=None, invoke_without_command=False, no_args_is_help=None, subcommand_metavar=None, chain=False, result_callback=None, **attrs)[source]

Bases: click.core.MultiCommand

Main cli class.

get_command(ctx, name)[source]

Return the requested command.

list_commands(ctx)[source]

Return the available commands (other than the indicators).

xclim.cli._create_command(indicator_name)[source]

Generate a Click.Command from an xclim Indicator.

xclim.cli._format_dict(data, formatter, key_fg='blue', spaces=2)[source]
xclim.cli._get_indicator(indicator_name)[source]
xclim.cli._get_input(ctx)[source]

Return the input dataset stored in the given context.

If the dataset is not open, opens it with open_dataset if a single path was given, or with open_mfdataset if a tuple or glob path was given.

xclim.cli._get_output(ctx)[source]

Return the output dataset stored in the given context.

If the output dataset doesn’t exist, create it.

xclim.cli._process_indicator(indicator, ctx, **params)[source]

Add given climate indicator to the output dataset from variables in the input dataset.

Computation is not triggered here if dask is enabled.

xclim.cli.write_file(ctx, *args, **kwargs)[source]

Write the output dataset to file.