Skip to content

Commit

Permalink
prob docstrings (#210)
Browse files Browse the repository at this point in the history
* add probabilistic dostrings
  • Loading branch information
aaronspring committed Oct 17, 2020
1 parent f628918 commit 4d271e8
Showing 1 changed file with 115 additions and 5 deletions.
120 changes: 115 additions & 5 deletions xskillscore/core/probabilistic.py
Original file line number Diff line number Diff line change
Expand Up @@ -55,6 +55,23 @@ def crps_gaussian(observations, mu, sig, dim=None, weights=None, keep_attrs=Fals
-------
xarray.Dataset or xarray.DataArray reduced by dimension dim
Examples
--------
>>> observations = xr.DataArray(np.random.normal(size=(3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3))])
>>> forecasts = xr.DataArray(np.random.normal(size=(3,3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3)),
... ('member', np.arange(3))])
>>> mu = forecasts.mean('member')
>>> sig = forecasts.std('member')
>>> crps_gaussian(observations, mu, sig, dim='x')
<xarray.DataArray (y: 3)>
array([1.0349773 , 0.36521376, 0.39017126])
Coordinates:
* y (y) int64 0 1 2
See Also
--------
properscoring.crps_gaussian
Expand Down Expand Up @@ -85,7 +102,7 @@ def crps_gaussian(observations, mu, sig, dim=None, weights=None, keep_attrs=Fals


def crps_quadrature(
x,
observations,
cdf_or_dist,
xmin=None,
xmax=None,
Expand All @@ -99,7 +116,7 @@ def crps_quadrature(
Parameters
----------
x : xarray.Dataset or xarray.DataArray
observations : xarray.Dataset or xarray.DataArray
Observations associated with the forecast distribution ``cdf_or_dist``.
cdf_or_dist : callable or scipy.stats.distribution
Function which returns the cumulative density of the forecast
Expand All @@ -121,13 +138,25 @@ def crps_quadrature(
-------
xarray.Dataset or xarray.DataArray
Examples
--------
>>> observations = xr.DataArray(np.random.normal(size=(3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3))])
>>> from scipy.stats import norm
>>> crps_quadrature(observations, norm)
<xarray.DataArray (y: 3)>
array([0.80280921, 0.31818197, 0.32364912])
Coordinates:
* y (y) int64 0 1 2
See Also
--------
properscoring.crps_quadrature
"""
res = xr.apply_ufunc(
properscoring.crps_quadrature,
x,
observations,
cdf_or_dist,
xmin,
xmax,
Expand All @@ -153,7 +182,7 @@ def crps_ensemble(
weights=None,
keep_attrs=False,
):
"""Continuous Ranked Probability Score with the ensemble distribution
"""Continuous Ranked Probability Score with the ensemble distribution.
Parameters
----------
Expand Down Expand Up @@ -187,6 +216,21 @@ def crps_ensemble(
-------
xarray.Dataset or xarray.DataArray
Examples
--------
>>> observations = xr.DataArray(np.random.normal(size=(3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3))])
>>> forecasts = xr.DataArray(np.random.normal(size=(3,3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3)),
... ('member', np.arange(3))])
>>> crps_ensemble(observations, forecasts, dim='x')
<xarray.DataArray (y: 3)>
array([1.04497153, 0.48997746, 0.47994095])
Coordinates:
* y (y) int64 0 1 2
See Also
--------
properscoring.crps_ensemble
Expand Down Expand Up @@ -237,6 +281,21 @@ def brier_score(observations, forecasts, dim=None, weights=None, keep_attrs=Fals
-------
xarray.Dataset or xarray.DataArray
Examples
--------
>>> observations = xr.DataArray(np.random.normal(size=(3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3))])
>>> forecasts = xr.DataArray(np.random.normal(size=(3,3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3)),
... ('member', np.arange(3))])
>>> brier_score(observations > .5, (forecasts > .5).mean('member'))
<xarray.DataArray (y: 3)>
array([0.51851852, 0.14814815, 0.37037037])
Coordinates:
* y (y) int64 0 1 2
See Also
--------
properscoring.brier_score
Expand Down Expand Up @@ -309,6 +368,25 @@ def threshold_brier_score(
observations. Otherwise, it will have an additional final dimension
corresponding to the threshold levels. Not implemented yet.)
Examples
--------
>>> observations = xr.DataArray(np.random.normal(size=(3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3))])
>>> forecasts = xr.DataArray(np.random.normal(size=(3,3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3)),
... ('member', np.arange(3))])
>>> threshold = [.2, .5, .8]
>>> threshold_brier_score(observations, forecasts, threshold)
<xarray.DataArray (y: 3, threshold: 3)>
array([[0.44444444, 0.51851852, 0.48148148],
[0.18518519, 0.14814815, 0.03703704],
[0.18518519, 0.37037037, 0.03703704]])
Coordinates:
* y (y) int64 0 1 2
* threshold (threshold) int64 1 2 3
See Also
--------
properscoring.threshold_brier_score
Expand Down Expand Up @@ -400,6 +478,22 @@ def rps(
xarray.Dataset or xarray.DataArray:
ranked probability score
Examples
--------
>>> observations = xr.DataArray(np.random.normal(size=(3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3))])
>>> forecasts = xr.DataArray(np.random.normal(size=(3,3,3)),
... coords=[('x', np.arange(3)),
... ('y', np.arange(3)),
... ('member', np.arange(3))])
>>> category_edges = np.array([.2, .5, .8])
>>> rps(observations, forecasts, category_edges)
<xarray.DataArray 'histogram_category' (y: 3)>
array([1. , 1. , 0.33333333])
Coordinates:
* y (y) int64 0 1 2
References
----------
https://www.cawcr.gov.au/projects/verification/verif_web_page.html#RPS
Expand Down Expand Up @@ -552,6 +646,14 @@ def discrimination(
>>> disc = discrimination(observed_event,
... forecast_event_likelihood,
... dim=['x','y'])
<xarray.DataArray (event: 2, forecast_probability: 5)>
array([[0., 1., 0., 0., 0.],
[0., 1., 0., 0., 0.]])
Coordinates:
* forecast_probability (forecast_probability) float64 0.1 0.3 0.5 0.7 0.9
* event (event) bool True False
References
----------
Expand Down Expand Up @@ -633,8 +735,16 @@ def reliability(
... coords=[('x', np.arange(3)),
... ('y', np.arange(3))])
>>> rel = reliability(observations > 0.1,
... (forecasts > 0.1).mean('ensemble'),
... (forecasts > 0.1).mean('member'),
... dim='x')
<xarray.DataArray (y: 3, forecast_probability: 5)>
array([[ nan, 1. , nan, nan, nan],
[0. , 0.5 , nan, nan, nan],
[ nan, 0.33333333, nan, nan, nan]])
Coordinates:
* y (y) int64 0 1 2
* forecast_probability (forecast_probability) float64 0.1 0.3 0.5 0.7 0.9
samples (y, forecast_probability) float64 0.0 3.0 ... 0.0 0.0
Notes
-----
Expand Down

0 comments on commit 4d271e8

Please sign in to comment.