Coverage for src/scores/probability/crps_impl.py: 100%
157 statements
« prev ^ index » next coverage.py v7.3.2, created at 2024-02-28 12:51 +1100
« prev ^ index » next coverage.py v7.3.2, created at 2024-02-28 12:51 +1100
1"""
2This module supports the implementation of the CRPS scoring function, drawing from additional functions.
3The two primary methods, `crps_cdf` and `crps_for_ensemble` are imported into
4the probability module to be part of the probability API.
5"""
6from collections.abc import Iterable
7from typing import Literal, Optional, Sequence
9import numpy as np
10import pandas as pd
11import xarray as xr
13import scores.utils
14from scores.probability.checks import coords_increasing
15from scores.probability.functions import (
16 add_thresholds,
17 cdf_envelope,
18 decreasing_cdfs,
19 integrate_square_piecewise_linear,
20 observed_cdf,
21 propagate_nan,
22)
25# pylint: disable=too-many-arguments
26# pylint: disable=too-many-branches
27def check_crps_cdf_inputs(
28 fcst,
29 obs,
30 threshold_dim,
31 threshold_weight,
32 fcst_fill_method,
33 threshold_weight_fill_method,
34 integration_method,
35 dims,
36):
37 """Checks that `crps_cdf` inputs are valid."""
38 if threshold_dim not in fcst.dims:
39 raise ValueError(f"'{threshold_dim}' is not a dimension of `fcst`")
41 if threshold_weight is not None and threshold_dim not in threshold_weight.dims:
42 raise ValueError(f"'{threshold_dim}' is not a dimension of `threshold_weight`")
44 if threshold_dim in obs.dims:
45 raise ValueError(f"'{threshold_dim}' is a dimension of `obs`")
47 if not set(obs.dims).issubset(fcst.dims):
48 raise ValueError("Dimensions of `obs` must be a subset of dimensions of `fcst`")
50 if threshold_weight is not None and not set(threshold_weight.dims).issubset(fcst.dims):
51 raise ValueError("Dimensions of `threshold_weight` must be a subset of dimensions of `fcst`")
53 if dims is not None and not set(dims).issubset(fcst.dims):
54 raise ValueError("`dims` must be a subset of `fcst` dimensions") # pragma: no cover
56 if fcst_fill_method not in ["linear", "step", "forward", "backward"]:
57 raise ValueError("`fcst_fill_method` must be 'linear', 'step', 'forward' or 'backward'")
59 if threshold_weight is not None and threshold_weight_fill_method not in [
60 "linear",
61 "step",
62 "forward",
63 "backward",
64 ]:
65 msg = "`threshold_weight_fill_method` must be 'linear', 'step', 'forward' or "
66 msg += "'backward' if `threshold_weight` is supplied"
67 raise ValueError(msg)
69 if integration_method not in ["exact", "trapz"]:
70 raise ValueError("`integration_method` must be 'exact' or 'trapz'")
72 if len(fcst[threshold_dim]) < 2:
73 raise ValueError("`threshold_dim` in `fcst` must have at least 2 values to calculate CRPS")
75 if not coords_increasing(fcst, threshold_dim):
76 raise ValueError("`threshold_dim` coordinates in `fcst` must be increasing")
78 if threshold_weight is not None and not coords_increasing(threshold_weight, threshold_dim):
79 raise ValueError("`threshold_dim` coordinates in `threshold_weight` must be increasing")
81 if threshold_weight is not None and (threshold_weight < 0).any():
82 raise ValueError("`threshold_weight` has negative values")
85def crps_cdf_reformat_inputs(
86 fcst: xr.DataArray,
87 obs: xr.DataArray,
88 threshold_dim: str,
89 threshold_weight: Optional[xr.DataArray],
90 additional_thresholds: Optional[Iterable[float]],
91 fcst_fill_method: Literal["linear", "step", "forward", "backward"],
92 threshold_weight_fill_method: Literal["linear", "step", "forward", "backward"],
93) -> tuple[xr.DataArray, xr.DataArray, xr.DataArray]:
94 """
95 Takes `fcst`, `obs` and `threshold_weight` inputs from `crps_cdf` and reformats them
96 for use in `crps_cdf`; that is, `cdf_fcst`, `cdf_obs` and `threshold_weight` are
97 returned with matching dimensions and coordinates.
99 `additional_thresholds` contains additional thresholds to use.
100 Useful if wanting to increase density of thresholds for `crps_cdf_trapz`.
102 NaNs are filled, unless there are too few non-NaNs for filling.
103 """
104 # will use all thresholds from fcst, obs and (if applicable) weight
106 fcst_thresholds = fcst[threshold_dim].values
107 obs_thresholds = pd.unique(obs.values.flatten())
109 weight_thresholds = []
110 if threshold_weight is not None:
111 weight_thresholds = threshold_weight[threshold_dim].values
113 if additional_thresholds is None:
114 additional_thresholds = []
116 thresholds = np.concatenate((weight_thresholds, fcst_thresholds, obs_thresholds, additional_thresholds))
117 thresholds = np.sort(pd.unique(thresholds))
118 thresholds = thresholds[~np.isnan(thresholds)]
120 # get obs in cdf form with correct thresholds
121 obs_cdf = observed_cdf(
122 obs,
123 threshold_dim,
124 threshold_values=thresholds,
125 include_obs_in_thresholds=False, # thresholds already has rounded obs values
126 precision=0,
127 )
129 # get fcst with correct thresholds
130 fcst = add_thresholds(fcst, threshold_dim, thresholds, fcst_fill_method)
132 # get weight with correct thresholds
133 if threshold_weight is None:
134 # the weight is 1
135 weight_cdf = xr.full_like(fcst, 1.0)
136 else:
137 weight_cdf = add_thresholds(threshold_weight, threshold_dim, thresholds, threshold_weight_fill_method)
139 fcst_cdf, obs_cdf, weight_cdf = xr.broadcast(fcst, obs_cdf, weight_cdf)
141 return fcst_cdf, obs_cdf, weight_cdf
144# pylint: disable=too-many-locals
145def crps_cdf(
146 fcst: xr.DataArray,
147 obs: xr.DataArray,
148 threshold_dim: str = "threshold",
149 threshold_weight: Optional[xr.DataArray] = None,
150 additional_thresholds: Optional[Iterable[float]] = None,
151 propagate_nans: bool = True,
152 fcst_fill_method: Literal["linear", "step", "forward", "backward"] = "linear",
153 threshold_weight_fill_method: Literal["linear", "step", "forward", "backward"] = "forward",
154 integration_method: Literal["exact", "trapz"] = "exact",
155 reduce_dims: Optional[Iterable[str]] = None,
156 preserve_dims: Optional[Iterable[str]] = None,
157 weights=None,
158 include_components=False,
159):
160 """Calculates the CRPS probabilistic metric given CDF input.
162 Calculates the continuous ranked probability score (CRPS), or the mean CRPS over
163 specified dimensions, given forecasts in the form of predictive cumulative
164 distribution functions (CDFs). Can also calculate threshold-weighted versions of the
165 CRPS by supplying a `threshold_weight`.
167 Predictive CDFs here are described by an indexed set of values rather than by
168 closed forumulae. As a result, the resolution or number of points at which the CDF
169 is realised has an impact on the calculation of areas under and over the curve to
170 obtain the CRPS result. The term 'threshold' is used to describe the dimension which
171 is used as an index for predictive CDF values. Various techniques are used to
172 interpolate CDF values between indexed thresholds.
174 Given
175 - a predictive CDF `fcst` indexed at thresholds by variable x,
176 - an observation in CDF form `obs_cdf` (i.e., obs_cdf(x) = 0 if x < obs and 1 if x >= obs),
177 - a `threshold_weight` array indexed by variable x,
179 The threshold-weighted CRPS is given by:
180 `CRPS = integral(threshold_weight(x) * (fcst(x) - obs_cdf(x))**2)`, over all thresholds x.
181 The usual CRPS is the threshold-weighted CRPS with `threshold_weight(x) = 1` for all x.
182 This can be decomposed into an over-forecast penalty
183 `integral(threshold_weight(x) * (fcst(x) - obs_cdf(x))**2)`, over all thresholds x where x >= obs
184 and an under-forecast penalty
185 `integral(threshold_weight(x) * (fcst(x) - obs_cdf(x))**2)`, over all thresholds x where x <= obs.
187 Note that the function `crps_cdf` is designed so that the `obs` argument contains
188 actual observed values. `crps_cdf` will convert `obs` into CDF form in order to
189 calculate the CRPS.
191 To calculate CRPS, integration is applied over the set of thresholds x taken from
192 - `fcst[threshold_dim].values`,
193 - `obs.values`.
194 - `threshold_weight[threshold_dim].values` if applicable.
195 - `additional_thresholds` if applicable.
196 With NaN values excluded. There are two methods of integration:
197 - "exact" gives the exact integral under that assumption that that `fcst` is
198 continuous and piecewise linear between its specified values, and that
199 `threshold_weight` (if supplied) is piecewise constant and right-continuous
200 between its specified values.
201 - "trapz" simply uses a trapezoidal rule using the specified values, and so is
202 an approximation of the CRPS. To get an accurate approximation, the density
203 of threshold values can be increased by supplying `additional_thresholds`.
205 Both methods of calculating CRPS may require adding additional values to the
206 `threshold_dim` dimension in `fcst` and (if supplied) `threshold_weight`.
207 `fcst_fill_method` and `weight_fill_method` specify how `fcst` and `weight` are to
208 be filled at these additional points.
210 The function `crps_cdf` calculates the CRPS using forecast CDF values 'as is'.
211 No checks are performed to ensure that CDF values in `fcst` are nondecreasing along
212 `threshold_dim`. Checks are conducted on `fcst` and `threshold_weight` (if applicable)
213 to ensure that coordinates are increasing along `threshold_dim`.
215 Args:
216 fcst: array of forecast CDFs, with the threshold dimension given by `threshold_dim`.
217 obs: array of observations, not in CDF form.
218 threshold_dim: name of the dimension in `fcst` that indexes the thresholds.
219 threshold_weight: weight to be applied along `threshold_dim` to calculat
220 threshold-weighted CRPS. Must contain `threshold_dim` as a dimension, and may
221 also include other dimensions from `fcst`. If `weight=None`, a weight of 1
222 is applied everywhere, which gives the usual CRPS.
223 additional_thresholds: additional thresholds values to add to `fcst` and (if
224 applicable) `threshold_weight` prior to calculating CRPS.
225 propagate_nans: If `True`, propagate NaN values along `threshold_dim` in `fcst`
226 and `threshold_weight` prior to calculating CRPS. This will result in CRPS
227 being NaN for these cases. If `False`, NaN values in `fcst` and `weight` will
228 be replaced, wherever possible, with non-NaN values using the fill method
229 specified by `fcst_fill_method` and `threshold_weight_fill_method`.
230 fcst_fill_method: how to fill values in `fcst` when NaNs have been introduced
231 (by including additional thresholds) or are specified to be removed (by
232 setting `propagate_nans=False`). Select one of:
233 - "linear": use linear interpolation, then replace any leading or
234 trailing NaNs using linear extrapolation. Afterwards, all values are
235 clipped to the closed interval [0, 1].
236 - "step": apply forward filling, then replace any leading NaNs with 0.
237 - "forward": first apply forward filling, then remove any leading NaNs by
238 back filling.
239 - "backward": first apply back filling, then remove any trailing NaNs by
240 forward filling.
241 In most cases, "linear" is likely the appropriate choice.
242 threshold_weight_fill_method: how to fill values in `threshold_weight` when NaNs
243 have been introduced (by including additional thresholds) or are specified
244 to be removed (by setting `propagate_nans=False`). Select one of "linear",
245 "step", "forward" or "backward". If the weight function is continuous,
246 "linear" is probably the best choice. If it is an increasing step function,
247 "forward" may be best.
248 integration_method (str): one of "exact" or "trapz".
249 preserve_dims (Tuple[str]): dimensions to preserve in the output. All other dimensions are collapsed
250 by taking the mean.
251 reduce_dims (Tuple[str]): dimensions to reduce in the output by taking the mean. All other dimensions are
252 preserved.
253 weights: Optionally provide an array for weighted averaging (e.g. by area, by latitude,
254 by population, custom)
255 include_components (bool): if True, include the under and over forecast components of
256 the score in the returned dataset.
258 Returns:
259 xr.Dataset: The following are the produced Dataset variables:
260 - "total" the total CRPS.
261 - "underforecast_penalty": the under-forecast penalty contribution of the CRPS.
262 - "overforecast_penalty": the over-forecast penalty contribution of the CRPS.
264 Raises:
265 ValueError: if `threshold_dim` is not a dimension of `fcst`.
266 ValueError: if `threshold_dim` is not a dimension of `threshold_weight` when
267 `threshold_weight` is not `None`.
268 ValueError: if `threshold_dim` is a dimension of `obs`.
269 ValueError: if dimensions of `obs` are not also dimensions of `fcst`.
270 ValueError: if dimensions of `threshold_weight` are not also dimensions of `fcst`
271 when `threshold_weight` is not `None`.
272 ValueError: if `dims` is not a subset of dimensions of `fcst`.
273 ValueError: if `fcst_fill_method` is not one of 'linear', 'step', 'forward' or
274 'backward'.
275 ValueError: if `weight_fill_method` is not one of 'linear', 'step', 'forward' or
276 'backward'.
277 ValueError: if `fcst[threshold_dim]` has less than 2 values.
278 ValueError: if coordinates in `fcst[threshold_dim]` are not increasing.
279 ValueError: if `threshold_weight` is not `None` and coordinates in
280 `threshold_weight[threshold_dim]` are not increasing.
281 ValueError: if `threshold_weight` has negative values.
283 See also:
284 - `scores.probability.crps_cdf_brier_decomposition`
285 - `scores.probability.crps_for_ensemble`
287 References:
288 - Matheson, J. E., and R. L. Winkler, 1976: Scoring rules for continuous probability distributions.
289 Manage. Sci.,22, 1087–1095.
290 - Gneiting, T., & Ranjan, R. (2011). Comparing Density Forecasts Using Threshold- and
291 Quantile-Weighted Scoring Rules.
292 Journal of Business & Economic Statistics, 29(3), 411–422. http://www.jstor.org/stable/23243806
293 """
295 dims = scores.utils.gather_dimensions(
296 fcst.dims,
297 obs.dims,
298 reduce_dims=reduce_dims,
299 preserve_dims=preserve_dims,
300 )
302 check_crps_cdf_inputs(
303 fcst,
304 obs,
305 threshold_dim,
306 threshold_weight,
307 fcst_fill_method,
308 threshold_weight_fill_method,
309 integration_method,
310 dims,
311 )
313 if propagate_nans:
314 fcst = propagate_nan(fcst, threshold_dim)
315 if threshold_weight is not None:
316 threshold_weight = propagate_nan(threshold_weight, threshold_dim)
318 fcst, obs, threshold_weight = crps_cdf_reformat_inputs(
319 fcst,
320 obs,
321 threshold_dim,
322 threshold_weight,
323 additional_thresholds,
324 fcst_fill_method,
325 threshold_weight_fill_method,
326 )
328 if integration_method == "exact":
329 result = crps_cdf_exact(
330 fcst,
331 obs,
332 threshold_weight,
333 threshold_dim,
334 include_components=include_components,
335 )
337 if integration_method == "trapz":
338 result = crps_cdf_trapz(
339 fcst,
340 obs,
341 threshold_weight,
342 threshold_dim,
343 include_components=include_components,
344 )
346 weighted = scores.functions.apply_weights(result, weights)
348 dims.remove(threshold_dim) # type: ignore
350 result = weighted.mean(dim=dims)
352 return result
355def crps_cdf_exact(
356 cdf_fcst: xr.DataArray,
357 cdf_obs: xr.DataArray,
358 threshold_weight: xr.DataArray,
359 threshold_dim: str,
360 include_components=False,
361) -> xr.Dataset:
362 """
363 Calculates exact value of CRPS assuming that:
364 - the forecast CDF is continuous piecewise linear, with join points given by
365 values in `cdf_fcst`,
366 - the observation CDF is right continuous with values in {0,1} given by `cdf_obs`,
367 - the threshold weight function is right continuous with values in {0,1} given
368 by `threshold_weight`.
370 If these assumptions do not hold, it might be best to use `crps_approximate`, with a
371 sufficiently high resolution along `threshold_dim`.
373 This function assumes that `cdf_fcst`, `cdf_obs`, `threshold_weight` have same shape.
374 Also assumes that values along the `threshold_dim` dimension are increasing.
376 Returns:
377 (xr.Dataset): Dataset with `threshold_dim` collapsed containing DataArrays with
378 CRPS and its decomposition, labelled "total", "underforecast_penalty" and
379 "overforecast_penalty". NaN is returned if there is a NaN in the corresponding
380 `cdf_fcst`, `cdf_obs` or `threshold_weight`.
381 """
383 # identify where input arrays have no NaN, collapsing `threshold_dim`
384 inputs_without_nan = (
385 ~np.isnan(cdf_fcst).any(threshold_dim)
386 & ~np.isnan(cdf_obs).any(threshold_dim)
387 & ~np.isnan(threshold_weight).any(threshold_dim)
388 )
390 # thresholds in the closure of the interval (i.e. including endpoints) where
391 # weight is 1
392 interval_where_weight_one = (threshold_weight == 1) | ((threshold_weight.shift(**{threshold_dim: 1})) == 1)
394 # thresholds in the closure of the interval where cdf_obs is 1
395 interval_where_obs_one = cdf_obs == 1
397 # thresholds in the closure of the interval where obs cdf is 0
398 interval_where_obs_zero = (cdf_obs == 0) | ((cdf_obs.shift(**{threshold_dim: 1})) == 0)
400 # over-forecast penalty contribution to CRPS: integral(w(x) * (F(x) - 1)^2) where x >= obs
401 over = (cdf_fcst - 1).where(interval_where_obs_one).where(interval_where_weight_one)
402 over = integrate_square_piecewise_linear(over, threshold_dim)
403 # If zero penalty, could be NaN. Replace with 0, then NaN using inputs_without_nan
404 over = over.where(~np.isnan(over), 0).where(inputs_without_nan)
406 # under-forecast penalty contribution to CRPS: integral(w(x) * F(x)^2) where x < obs
407 under = (cdf_fcst).where(interval_where_obs_zero).where(interval_where_weight_one)
408 under = integrate_square_piecewise_linear(under, threshold_dim)
409 # If zero penalty, could be NaN. Replace with 0, then nan using inputs_without_nan
410 under = under.where(~np.isnan(under), 0).where(inputs_without_nan)
412 total = over + under
413 result = total.to_dataset(name="total")
415 if include_components:
416 result = xr.merge(
417 [
418 total.rename("total"),
419 under.rename("underforecast_penalty"),
420 over.rename("overforecast_penalty"),
421 ]
422 )
424 return result
427def crps_cdf_brier_decomposition(
428 fcst: xr.DataArray,
429 obs: xr.DataArray,
430 threshold_dim: str = "threshold",
431 additional_thresholds: Optional[Iterable[float]] = None,
432 fcst_fill_method: Literal["linear", "step", "forward", "backward"] = "linear",
433 reduce_dims: Optional[Iterable[str]] = None,
434 preserve_dims: Optional[Iterable[str]] = None,
435) -> xr.Dataset:
436 """
437 Given an array `fcst` of predictive CDF values indexed along `threshold_dim`, and
438 an array `obs` of observations, calculates the mean Brier score for each index along
439 `threshold_dim`. Since the mean CRPS is the integral of the mean Brier score over
440 all thresholds, this gives a threshold decomposition of the mean CRPS.
442 If any there are any NaNs along the threshold dimension of `fcst`, then NaNs are
443 propagated along this dimension prior to calculating the decomposition. If
444 propagating NaNs is not desired, the user may first fill NaNs in `fcst` using
445 `scores.probability.functions.fill_cdf`.
447 Args:
448 fcst (xr.DataArray): DataArray of CDF values with threshold dimension `threshold_dim`.
449 obs (xr.DataArray): DataArray of observations, not in CDF form.
450 threshold_dim (str): name of the threshold dimension in `fcst`.
451 additional_thresholds (Optional[Iterable[float]]): additional thresholds
452 at which to calculate the mean Brier score.
453 fcst_fill_method (Literal["linear", "step", "forward", "backward"]): How to fill NaN
454 values in `fcst` that arise from new user-supplied thresholds or thresholds derived
455 from observations.
456 - "linear": use linear interpolation, and if needed also extrapolate linearly.
457 Clip to 0 and 1. Needs at least two non-NaN values for interpolation,
458 so returns NaNs where this condition fails.
459 - "step": use forward filling then set remaining leading NaNs to 0.
460 Produces a step function CDF (i.e. piecewise constant).
461 - "forward": use forward filling then fill any remaining leading NaNs with
462 backward filling.
463 - "backward": use backward filling then fill any remaining trailing NaNs with
464 forward filling.
465 dims: dimensions to preserve in the output. The dimension `threshold_dim` is always
466 preserved, even if not specified here.
468 Returns:
469 An xarray Dataset with data_vars:
470 - "total_penalty": the mean Brier score for each threshold.
471 - "underforecast_penalty": the mean of the underforecast penalties for the Brier score. For a particular
472 forecast case, this component equals 0 if the event didn't occur and the Brier score if it did.
473 - "overforecast_penalty": the mean of the overforecast penalties for the Brier score. For a particular
474 forecast case, this component equals 0 if the event did occur and the Brier score if it did not.
476 Raises:
477 ValueError: if `threshold_dim` is not a dimension of `fcst`.
478 ValueError: if `threshold_dim` is a dimension of `obs`.
479 ValueError: if dimensions of `obs` are not also among the dimensions of `fcst`.
480 ValueError: if dimensions in `dims` is not among the dimensions of `fcst`.
481 ValueError: if `fcst_fill_method` is not one of 'linear', 'step', 'forward' or 'backward'.
482 ValueError: if coordinates in `fcst[threshold_dim]` are not increasing.
483 """
484 dims = scores.utils.gather_dimensions(
485 fcst.dims,
486 obs.dims,
487 reduce_dims=reduce_dims,
488 preserve_dims=preserve_dims,
489 )
491 check_crps_cdf_brier_inputs(fcst, obs, threshold_dim, fcst_fill_method, dims)
493 fcst = propagate_nan(fcst, threshold_dim)
495 fcst, obs, _ = crps_cdf_reformat_inputs(
496 fcst,
497 obs,
498 threshold_dim,
499 None,
500 additional_thresholds,
501 fcst_fill_method,
502 "forward",
503 )
505 dims.remove(threshold_dim) # type: ignore
507 # brier score for each forecast case
508 bscore = (fcst - obs) ** 2
509 not_nan = ~np.isnan(bscore)
511 # `obs` here is the empirical CDF of the observation
512 # when `obs == 1` the observation was lower than the threshold considered
513 # and hence the Brier score at this threshold penalises an over-forecast
514 over = bscore.where(np.isclose(obs, 1), 0).where(not_nan).mean(dims)
515 under = bscore.where(np.isclose(obs, 0), 0).where(not_nan).mean(dims)
516 total = over + under
518 result = xr.merge(
519 [
520 total.rename("total_penalty"),
521 under.rename("underforecast_penalty"),
522 over.rename("overforecast_penalty"),
523 ]
524 )
526 return result
529def check_crps_cdf_brier_inputs(fcst, obs, threshold_dim, fcst_fill_method, dims):
530 """Checks inputs to 'crps_cdf_brier_decomposition'."""
532 if threshold_dim not in fcst.dims:
533 raise ValueError(f"'{threshold_dim}' is not a dimension of `fcst`")
535 if threshold_dim in obs.dims:
536 raise ValueError(f"'{threshold_dim}' is a dimension of `obs`")
538 if not set(obs.dims).issubset(fcst.dims):
539 raise ValueError("Dimensions of `obs` must be a subset of dimensions of `fcst`")
541 if dims is not None and not set(dims).issubset(fcst.dims):
542 raise ValueError("`dims` must be a subset of `fcst` dimensions") # pragma: no cover
544 if fcst_fill_method not in ["linear", "step", "forward", "backward"]:
545 raise ValueError("`fcst_fill_method` must be 'linear', 'step', 'forward' or 'backward'")
547 if not coords_increasing(fcst, threshold_dim):
548 raise ValueError("`threshold_dim` coordinates in `fcst` must be increasing")
551def adjust_fcst_for_crps(
552 fcst: xr.DataArray,
553 threshold_dim: str,
554 obs: xr.DataArray,
555 decreasing_tolerance: float = 0,
556 additional_thresholds: Optional[Iterable[float]] = None,
557 fcst_fill_method: Literal["linear", "step", "forward", "backward"] = "linear",
558 integration_method: Literal["exact", "trapz"] = "exact",
559) -> xr.DataArray:
560 """
561 This function takes a forecast cumulative distribution functions (CDF) `fcst`.
562 If `fcst` is not decreasing outside of specified tolerance, it returns `fcst`.
563 Otherwise, the CDF envelope for `fcst` is computed, and the CDF from among
564 - `fcst`,
565 - the upper envelope, and
566 - the lower envelope
567 that has the higher (i.e. worse) CRPS is returned. In the event of a tie,
568 preference is given in the order `fcst` then upper.
569 See `scores.probability.functions.cdf_envelope` for details about the CDF envelope.
571 The use case for this is when, either due to rounding or poor forecast process, the
572 forecast CDF `fcst` fails to be nondecreasing. Rather than simply replacing `fcst`
573 with NaN, `adjust_fcst_for_crps` returns a CDF for which CRPS can be calculated, but
574 possibly with a predictive performance cost as measured by CRPS.
576 Whether a CDF is decreasing outside specified tolerance is determined as follows.
577 For each CDF in `fcst`, the sum of incremental decreases along the threshold dimension
578 is calculated. For example, if the CDF values are
579 [0, 0.4, 0.3, 0.9, 0.88, 1]
580 then the sum of incremental decreases is -0.12. This CDF decreases outside specified
581 tolerance if 0.12 > `decreasing_tolerance`.
583 The adjusted array of forecast CDFs is determined as follows:
584 - any NaN values in `fcst` are propagated along `threshold_dim` so that in each case
585 the entire CDF is NaN;
586 - any CDFs in `fcst` that are decreasing within specified tolerance are unchanged;
587 - any CDFs in `fcst` that are decreasing outside specified tolerance are replaced with
588 whichever of the upper or lower CDF envelope gives the highest CRPS, unless the original
589 values give a higher CRPS in which case original values are kept.
590 See `scores.probability.functions.cdf_envelope` for a description of the 'CDF envelope'.
591 If propagating NaNs is not desired, the user may first fill NaNs in `fcst` using
592 `scores.probability.functions.fill_cdf`.
593 The CRPS for each forecast case is calculated using `crps`, with a weight of 1.
595 Args:
596 fcst: DataArray of CDF values with threshold dimension `threshold_dim`.
597 threshold_dim: name of the threshold dimension in `fcst`.
598 obs: DataArray of observations.
599 decreasing_tolerance: nonnegative tolerance value.
600 additional_thresholds: optional additional thresholds passed on to `crps` when
601 calculating CRPS.
602 fcst_fill_method: `fcst` fill method passed on to `crps` when calculating CRPS.
603 integration_method: integration method passed on to `crps` when calculating CRPS.
605 Returns:
606 An xarray DataArray of possibly adjusted forecast CDFs, where adjustments are made
607 to penalise CDFs that decrease outside tolerance.
609 Raises:
610 ValueError: If `threshold_dim` is not a dimension of `fcst`.
611 ValueError: If `decreasing_tolerance` is negative.
612 """
613 if threshold_dim not in fcst.dims:
614 raise ValueError(f"'{threshold_dim}' is not a dimension of `fcst`")
616 if decreasing_tolerance < 0:
617 raise ValueError("`decreasing_tolerance` must be nonnegative")
619 fcst = propagate_nan(fcst, threshold_dim)
621 is_decreasing = decreasing_cdfs(fcst, threshold_dim, decreasing_tolerance)
623 if not is_decreasing.any():
624 return fcst
626 fcst_env = cdf_envelope(fcst, threshold_dim)
628 # dimensions to preserve when calculating CRPS
629 crps_dims = [x for x in fcst_env.dims if x != threshold_dim]
631 crps_fcst_env = crps_cdf(
632 fcst_env,
633 obs,
634 threshold_dim,
635 additional_thresholds=additional_thresholds,
636 fcst_fill_method=fcst_fill_method,
637 integration_method=integration_method,
638 preserve_dims=crps_dims,
639 )
641 fcst_type_to_use = crps_fcst_env["total"].idxmax("cdf_type")
643 fcst_to_return = fcst_env.where(fcst_env["cdf_type"] == fcst_type_to_use).max("cdf_type")
645 fcst_to_return = fcst_to_return.where(is_decreasing).combine_first(fcst)
647 return fcst_to_return
650def crps_cdf_trapz(
651 cdf_fcst: xr.DataArray,
652 cdf_obs: xr.DataArray,
653 threshold_weight: xr.DataArray,
654 threshold_dim: str,
655 include_components=False,
656) -> xr.Dataset:
657 """
658 Returns dataset with CRPS estimate and decomposition using a trapezoidal rule on
659 points from the integrand
660 threshold_weight(x) * (cdf_fcst(x) - cdf_obs(x)) ** 2
662 The error in the estimate improves as the distance between adjacent thresholds in
663 `threshold_dim` approaches zero.
665 NaN is returned if there is a NaN in the corresponding `cdf_fcst`, `cdf_obs` or
666 `threshold_weight`.
667 """
669 # identify where input arrays have no NaN, collapsing `threshold_dim`
670 inputs_without_nan = (
671 ~np.isnan(cdf_fcst).any(threshold_dim)
672 & ~np.isnan(cdf_obs).any(threshold_dim)
673 & ~np.isnan(threshold_weight).any(threshold_dim)
674 )
676 # total error measured by CRPS
677 total = (threshold_weight * (cdf_fcst - cdf_obs) ** 2).integrate(threshold_dim).where(inputs_without_nan)
679 over = (cdf_obs * threshold_weight * (cdf_fcst - cdf_obs) ** 2).integrate(threshold_dim).where(inputs_without_nan)
681 under = total - over
683 result = total.to_dataset(name="total")
685 if include_components:
686 result = xr.merge(
687 [
688 total.rename("total"),
689 under.rename("underforecast_penalty"),
690 over.rename("overforecast_penalty"),
691 ]
692 )
694 return result
697def crps_step_threshold_weight(
698 step_points: xr.DataArray,
699 threshold_dim: str,
700 threshold_values: Optional[Iterable[float]] = None,
701 steppoints_in_thresholds: bool = True,
702 steppoint_precision: float = 0,
703 weight_upper: bool = True,
704) -> xr.DataArray:
705 """Generates an array of weights based on DataArray step points.
707 Creates an array of threshold weights, which can be used to calculate
708 threshold-weighted CRPS, based on a step function. Applies a weight of 1 when
709 step_point >= threshold, and a weight of 0 otherwise. Zeros and ones in the output
710 weight function can be reversed by setting `weight_upper=False`.
712 Args:
713 step_points (xr.DataArray): points at which the weight function changes value from 0 to 1.
714 threshold_dim (str): name of the threshold dimension in the returned array weight function.
715 threshold_values (str): thresholds at which to calculate weights.
716 steppoints_in_thresholds (bool): include `step_points` among the `threshold_dim` values.
717 steppoint_precision (float): precision at which to round step_points prior to calculating the
718 weight function. Select 0 for no rounding.
719 weight_upper (bool): If true, returns a weight of 1 if step_point >= threshold, and a
720 weight of 0 otherwise. If false, returns a weight of 0 if step_point >= threshold,
721 and a weight of 1 otherwise.
723 Returns:
724 (xr.DataArray): Zeros and ones with the dimensions in `step_points`
725 and an additional `threshold_dim` dimension.
726 """
728 weight = observed_cdf(
729 step_points,
730 threshold_dim,
731 threshold_values=threshold_values,
732 include_obs_in_thresholds=steppoints_in_thresholds,
733 precision=steppoint_precision,
734 )
736 if not weight_upper:
737 weight = 1 - weight
739 return weight
742def crps_for_ensemble(
743 fcst: xr.DataArray,
744 obs: xr.DataArray,
745 ensemble_member_dim: str,
746 method: Literal["ecdf", "fair"] = "ecdf",
747 reduce_dims: Optional[Sequence[str]] = None,
748 preserve_dims: Optional[Sequence[str]] = None,
749 weights: xr.DataArray = None,
750) -> xr.DataArray:
751 """Calculates the CRPS probabilistic metric given ensemble input.
753 Calculates the continuous ranked probability score (CRPS) given an ensemble of forecasts.
754 An ensemble of forecasts can also be thought of as a random sample from the predictive
755 distribution.
757 Given an observation y, and ensemble member values {x_i} (for 1 <= i <= M), the CRPS is
758 calculated by the formula
759 CRPS({x_i}, y) = (1 / M) * sum(|x_i - y|) - (1 / 2 * K) * sum(|x_i - x_j|),
760 where the first sum is iterated over 1 <= i <= M and the second sum is iterated over
761 1 <= i <= M and 1 <= j <= M.
763 The value of the constant K in this formula depends on the method.
764 - If `method="ecdf"` then K = M ** 2. In this case the CRPS value returned is
765 the exact CRPS value for the emprical cumulation distribution function
766 constructed using the ensemble values.
767 - If `method="fair"` then K = M * (M - 1). In this case the CRPS value returned
768 is the approximated CRPS where the ensemble values can be interpreted as a
769 random sample from the underlying predictive distribution. This interpretation
770 stems from the formula CRPS(F, Y) = E|X - Y| - E|X - X'|/2, where X and X'
771 are independent samples of the predictive distribution F, Y is the observation
772 (possibly unknown) and E denotes the expectation. This choice of K gives an
773 unbiased estimate for the second expectation.
775 Args:
776 fcst: Forecast data. Must have a dimension `ensemble_member_dim`.
777 obs: Observation data.
778 ensemble_member_dim: the dimension that specifies the ensemble member or the sample
779 from the predictive distribution.
780 method: Either "ecdf" or "fair".
781 reduce_dims: Dimensions to reduce. Can be "all" to reduce all dimensions.
782 preserve_dims: Dimensions to preserve. Can be "all" to preserve all dimensions.
783 weights: Weights for calculating a weighted mean of individual scores.
785 Returns:
786 xarray object of (weighted mean) CRPS values.
788 Raises:
789 ValueError: when method is not one of "ecdf" or "fair".
791 See also:
792 `scores.probability.crps_cdf`
794 References:
795 - C. Ferro (2014), "Fair scores for ensemble forecasts", Q J R Meteorol Soc
796 140(683):1917-1923.
797 - T. Gneiting T and A. Raftery (2007), "Strictly proper scoring rules, prediction,
798 and estimation", J Am Stat Assoc, 102(477):359-37.
799 - M. Zamo and P. Naveau (2018), "Estimation of the Continuous Ranked Probability
800 Score with Limited Information and Applications to Ensemble Weather Forecasts",
801 Math Geosci 50:209-234, https://doi.org/10.1007/s11004-017-9709-7
802 """
803 if method not in ["ecdf", "fair"]:
804 raise ValueError("`method` must be one of 'ecdf' or 'fair'")
806 dims_for_mean = scores.utils.gather_dimensions2(fcst, obs, weights, reduce_dims, preserve_dims, ensemble_member_dim)
808 ensemble_member_dim1 = scores.utils.tmp_coord_name(fcst)
810 # calculate forecast spread contribution
811 fcst_copy = fcst.rename({ensemble_member_dim: ensemble_member_dim1})
813 fcst_spread_term = np.abs(fcst - fcst_copy).sum([ensemble_member_dim, ensemble_member_dim1])
814 ens_count = fcst.count(ensemble_member_dim)
815 if method == "ecdf":
816 fcst_spread_term = fcst_spread_term / (2 * ens_count**2)
817 if method == "fair":
818 fcst_spread_term = fcst_spread_term / (2 * ens_count * (ens_count - 1))
820 # calculate final CRPS for each forecast case
821 fcst_obs_term = np.abs(fcst - obs).mean(ensemble_member_dim)
822 result = fcst_obs_term - fcst_spread_term
824 # apply weights and take means across specified dims
825 result = scores.functions.apply_weights(result, weights).mean(dim=dims_for_mean)
827 return result