Coverage for src/scores/continuous/murphy_impl.py: 100%
88 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"""
2Murphy score
3"""
4from collections.abc import Sequence
5from typing import Literal, Optional, Union
7import numpy as np
8import xarray as xr
10from scores.processing import broadcast_and_match_nan
11from scores.typing import FlexibleArrayType, FlexibleDimensionTypes
12from scores.utils import gather_dimensions
14QUANTILE = "quantile"
15HUBER = "huber"
16EXPECTILE = "expectile"
17VALID_SCORING_FUNC_NAMES = [QUANTILE, HUBER, EXPECTILE]
18SCORING_FUNC_DOCSTRING_PARAMS = {fun.upper(): fun for fun in VALID_SCORING_FUNC_NAMES}
21def murphy_score(
22 fcst: xr.DataArray,
23 obs: xr.DataArray,
24 thetas: Union[Sequence[float], xr.DataArray],
25 functional: Literal["quantile", "huber", "expectile"],
26 alpha: float,
27 huber_a: Optional[float] = None,
28 decomposition: bool = False,
29 reduce_dims: FlexibleDimensionTypes = None,
30 preserve_dims: FlexibleDimensionTypes = None,
31) -> xr.Dataset:
32 """Returns the mean elementary score (Ehm et. al. 2016), also known as Murphy score,
33 evaluated at decision thresholds specified by thetas. Optionally returns a decomposition
34 of the score in terms of penalties for over- and under-forecasting.
36 Select `functional="quantile"` and `alpha=0.5` for the median functional.
37 Select `functional="expectile"` and `alpha=0.5` for the mean (i.e., expectation) functional.
39 Consider using `murphy_thetas` to generate thetas. If utilising dask, you may want
40 to store `thetas` as a netCDF on disk and pass it in as an xarray object. Otherwise,
41 very large objects may be created when `fcst`, `obs` and `thetas` are broadcast
42 together.
45 Args:
46 fcst: Forecast numerical values.
47 obs: Observed numerical values.
48 thetas: Theta thresholds.
49 functional: The type of elementary scoring function, one of {QUANTILE},
50 {HUBER}, or {EXPECTILE}.
51 alpha: Risk parameter (i.e. quantile or expectile level) for the functional. Must be between 0 and 1.
52 huber_a: Huber transition parameter, used for {HUBER} functional only.
53 Must be strictly greater than 0.
54 decomposition: True to return penalty values due to under- and over-fcst
55 as well as the total score, False to return total score only.
56 reduce_dims: Optionally specify which dimensions to reduce when
57 calculating the Murphy score. All other dimensions will be preserved. As a
58 special case, 'all' will allow all dimensions to be reduced. Only one
59 of `reduce_dims` and `preserve_dims` can be supplied. The default behaviour
60 if neither are supplied is to reduce all dims.
61 preserve_dims: Optionally specify which dimensions to preserve
62 when calculating the Murphy score. All other dimensions will be reduced.
63 As a special case, 'all' will allow all dimensions to be
64 preserved. In this case, the result will be in the same
65 shape/dimensionality as the forecast, and the errors will be
66 the FIRM score at each point (i.e. single-value comparison
67 against observed), and the forecast and observed dimensions
68 must match precisely. Only one of `reduce_dims` and `preserve_dims` can be
69 supplied. The default behaviour if neither are supplied is to reduce all dims.
71 Returns:
72 An xr.Dataset with dimensions based on the `preserve_dims` or `reduce_dims` arg
73 as well as a "theta" dimension with values matching `thetas` input.
75 If `decomposition` is False, the dataset's variables will contain 1 element
76 "total".
78 If `decomposition` is True, in addition to "total", it will have
79 "underforecast" and "overforecast" data_vars.
81 Raises:
82 ValueError: If `functional` is not one of the expected functions.
83 ValueError: If `alpha` is not strictly between 0 and 1.
84 ValueError: If `huber_a` is not strictly greater than 0.
86 References:
87 For mean elementary score definitions, see
89 - Theorem 1 of Ehm et. al. (2016) "Of Quantiles and Expectiles", J. Royal Stat. Soc. B,
90 78(3): 505–562. https://www.jstor.org/stable/24775351.
91 - Theorem 5.3 of Taggart (2022) "Point forecasting and forecast evaluation with
92 generalized Huber loss", Electron. J. Statist. 16(1): 201-231.
93 DOI: 10.1214/21-EJS1957
95 """
96 functional_lower = functional.lower()
97 _check_murphy_inputs(alpha=alpha, functional=functional_lower, huber_a=huber_a)
98 if isinstance(thetas, xr.DataArray):
99 theta1 = thetas
100 else:
101 theta1 = xr.DataArray(data=thetas, dims=["theta"], coords=dict(theta=thetas))
102 theta1, fcst1, obs1 = broadcast_and_match_nan(theta1, fcst, obs)
104 over, under = globals()[f"_{functional_lower}_elementary_score"](fcst1, obs1, theta1, alpha, huber_a=huber_a)
105 # Align dimensions, this is required in cases such as when the station numbers
106 # are not in the same order in `obs` and `fcst` to prevent an exception on the next
107 # line that combines the scores
108 over, under, fcst1, *_ = xr.align(over, under, fcst1)
109 score = over.combine_first(under).fillna(0).where(~np.isnan(fcst1), np.nan)
111 sources = [score]
112 names = ["total"]
113 if decomposition:
114 over = over.fillna(0).where(~np.isnan(fcst1), np.nan)
115 under = under.fillna(0).where(~np.isnan(fcst1), np.nan)
116 sources += [under, over]
117 names += ["underforecast", "overforecast"]
118 for source, name in zip(sources, names):
119 source.name = name
120 result = xr.merge(sources)
121 reduce_dims = gather_dimensions(fcst.dims, obs.dims, reduce_dims, preserve_dims)
122 result = result.mean(dim=reduce_dims)
123 return result
126def _quantile_elementary_score(fcst: FlexibleArrayType, obs: FlexibleArrayType, theta, alpha, **_):
127 """Return over and under forecast vs obs penalties relative to theta for {QUANTILE}."""
128 zero_array = fcst * 0.0
129 over = (zero_array + 1 - alpha).where((fcst > theta) & (obs <= theta))
130 under = (zero_array + alpha).where((fcst <= theta) & (obs > theta))
131 return over, under
134def _huber_elementary_score(fcst: FlexibleArrayType, obs: FlexibleArrayType, theta, alpha, *, huber_a):
135 """Return over and under forecast vs obs penalties relative to theta for {HUBER}."""
136 zero_array = fcst * 0.0
137 over = (1 - alpha) * np.minimum(theta - obs, zero_array + huber_a).where((fcst > theta) & (obs <= theta))
138 under = alpha * np.minimum(obs - theta, zero_array + huber_a).where((fcst <= theta) & (obs > theta))
139 return over, under
142def _expectile_elementary_score(fcst: FlexibleArrayType, obs: FlexibleArrayType, theta, alpha, **_):
143 """Return over and under forecast vs obs penalties relative to theta for {EXPECTILE}."""
144 over = (1 - alpha) * np.abs(obs - theta).where((fcst > theta) & (obs <= theta))
145 under = alpha * np.abs(obs - theta).where((fcst <= theta) & (obs > theta))
146 return over, under
149def _check_murphy_inputs(alpha=None, functional=None, huber_a=None, left_limit_delta=None):
150 """Raise ValueError if the arguments have unexpected values."""
151 if (alpha is not None) and not (0 < alpha < 1):
152 err = f"alpha (={alpha}) argument for Murphy scoring function should be strictly " "between 0 and 1."
153 raise ValueError(err)
154 if (functional is not None) and (functional not in VALID_SCORING_FUNC_NAMES):
155 err = (
156 f"Functional option '{functional}' for Murphy scoring function is "
157 f"unknown, should be one of {VALID_SCORING_FUNC_NAMES}."
158 )
159 raise ValueError(err)
160 if (functional == HUBER) and ((huber_a is None) or (huber_a <= 0)):
161 err = f"huber_a (={huber_a}) argument should be > 0 when functional='{HUBER}'."
162 raise ValueError(err)
163 if (left_limit_delta is not None) and (left_limit_delta < 0):
164 err = f"left_limit_delta (={left_limit_delta}) argument should be >= 0."
165 raise ValueError(err)
168def murphy_thetas(
169 forecasts: list[xr.DataArray],
170 obs: xr.DataArray,
171 functional: Literal["quantile", "huber", "expectile"],
172 huber_a: Optional[float] = None,
173 left_limit_delta: Optional[float] = None,
174) -> list[float]:
175 """Return the decision thresholds (theta values) at which to evaluate
176 elementary scores for the construction of Murphy diagrams.
178 This function may generate a very large number of theta thresholds if the forecast
179 and obs data arrays are large and their values have high precision, which may lead
180 to long computational times for Murphy diagrams. To reduce the number of thetas,
181 users can first round forecast and obs values to an appropriate resolution.
183 Args:
184 forecasts: Forecast values, one array per source.
185 obs: Observed values.
186 functional: The type of elementary scoring function, one of {QUANTILE},
187 {HUBER}, or {EXPECTILE}.
188 huber_a: Huber transition parameter, used for {HUBER} functional only.
189 Must be strictly greater than 0.
190 left_limit_delta: Approximation of the left hand limit, used for {HUBER}
191 and {EXPECTILE} functionals. Must be greater than or equal to 0.
192 None will be treated as 0. Ideally, left_limit_delta should be
193 small relative to the fcst and obs precision, and not greater than
194 that precision.
196 Returns:
197 List[float]: theta thresholds to be used to compute murphy scores.
199 Raises:
200 ValueError: If `functional` is not one of the expected functions.
201 ValueError: If `huber_a` is not strictly greater than 0.
202 ValueError: If `left_limit_delta` is not greater than or equal to 0.
204 Notes:
205 For theta values at which evaluate elementary scores, see
207 - Corollary 2 (p.521) of Ehm et. al. (2016) "Of Quantiles and Expectiles", J. Royal Stat. Soc. B,
208 78(3): 505–562. https://www.jstor.org/stable/24775351.
209 - Corollary 5.6 of Taggart (2022) "Point forecasting and forecast evaluation with
210 generalized Huber loss", Electron. J. Statist. 16(1): 201-231.
211 DOI: 10.1214/21-EJS1957
213 """
214 _check_murphy_inputs(functional=functional, huber_a=huber_a, left_limit_delta=left_limit_delta)
215 if (left_limit_delta is None) and (functional in ["huber", "expectile"]):
216 left_limit_delta = 0
218 func = globals()[f"_{functional}_thetas"]
219 result = func(
220 forecasts=forecasts,
221 obs=obs,
222 huber_a=huber_a,
223 left_limit_delta=left_limit_delta,
224 )
225 return result
228def _quantile_thetas(forecasts, obs, **_):
229 """Return thetas for {QUANTILE} elementary scoring function."""
230 ufcasts_and_uobs = np.unique(np.concatenate([*[x.values.flatten() for x in forecasts], obs.values.flatten()]))
231 result = ufcasts_and_uobs[~np.isnan(ufcasts_and_uobs)]
232 return list(result)
235def _huber_thetas(forecasts, obs, *, huber_a, left_limit_delta, **_):
236 """Return thetas for {HUBER} elementary scoring function."""
237 uobs = np.unique(obs)
238 uobs_minus_a = uobs - huber_a
239 uobs_plus_a = uobs + huber_a
240 ufcasts = np.unique(forecasts)
241 left_limit_points = ufcasts - left_limit_delta
242 ufcasts_and_uobs = np.unique(np.concatenate([ufcasts, left_limit_points, uobs, uobs_minus_a, uobs_plus_a]))
243 result = ufcasts_and_uobs[~np.isnan(ufcasts_and_uobs)]
244 return list(result)
247def _expectile_thetas(forecasts, obs, *, left_limit_delta, **_):
248 """Return thetas for {EXPECTILE} elementary scoring function."""
249 ufcasts = np.unique(forecasts)
250 left_limit_points = ufcasts - left_limit_delta
251 ufcasts_and_uobs = np.unique(np.concatenate([ufcasts, left_limit_points, obs.values.flatten()]))
252 result = ufcasts_and_uobs[~np.isnan(ufcasts_and_uobs)]
253 return list(result)