Skip to content

Efficiency Metrics

Model efficiency and performance measures for evaluating forecast skill.

Efficiency Metrics for Model Evaluation (Aero Protocol Compliant).

KGE(obs, mod, axis=None)

Kling-Gupta Efficiency (KGE).

Typical Use Cases

  • Quantifying the overall agreement between model and observations, combining correlation, bias, and variability.
  • Used in hydrology, meteorology, and environmental model evaluation.

Typical Values and Range

  • Range: -∞ to 1
  • 1: Perfect agreement between model and observations
  • 0: Moderate skill
  • Negative values: Poor skill

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values. mod : numpy.ndarray or xarray.DataArray Model or predicted values. axis : int, str, or iterable of such, optional Axis along which to compute KGE.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Kling-Gupta efficiency (unitless, -∞ to 1).

Examples

import numpy as np from monet_stats.correlation_metrics import KGE obs = np.array([1, 2, 3]) mod = np.array([1.1, 1.9, 3.2]) KGE(obs, mod) 0.8988771192996924

Source code in src/monet_stats/correlation_metrics.py
 921
 922
 923
 924
 925
 926
 927
 928
 929
 930
 931
 932
 933
 934
 935
 936
 937
 938
 939
 940
 941
 942
 943
 944
 945
 946
 947
 948
 949
 950
 951
 952
 953
 954
 955
 956
 957
 958
 959
 960
 961
 962
 963
 964
 965
 966
 967
 968
 969
 970
 971
 972
 973
 974
 975
 976
 977
 978
 979
 980
 981
 982
 983
 984
 985
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
def KGE(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Kling-Gupta Efficiency (KGE).

    Typical Use Cases
    -----------------
    - Quantifying the overall agreement between model and observations,
      combining correlation, bias, and variability.
    - Used in hydrology, meteorology, and environmental model evaluation.

    Typical Values and Range
    ------------------------
    - Range: -∞ to 1
    - 1: Perfect agreement between model and observations
    - 0: Moderate skill
    - Negative values: Poor skill

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values.
    mod : numpy.ndarray or xarray.DataArray
        Model or predicted values.
    axis : int, str, or iterable of such, optional
        Axis along which to compute KGE.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Kling-Gupta efficiency (unitless, -∞ to 1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import KGE
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([1.1, 1.9, 3.2])
    >>> KGE(obs, mod)
    0.8988771192996924
    """
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        # Handle axis vs dim
        if axis is not None:
            if isinstance(axis, int):
                dim = obs.dims[axis]
            elif isinstance(axis, (list, tuple)):
                dim = [obs.dims[d] if isinstance(d, int) else d for d in axis]
            else:
                dim = axis
        else:
            dim = obs.dims

        r = xr.corr(obs, mod, dim=dim)
        alpha = mod.std(dim=dim) / obs.std(dim=dim)
        beta = mod.mean(dim=dim) / obs.mean(dim=dim)
        result = 1.0 - ((r - 1.0) ** 2 + (alpha - 1.0) ** 2 + (beta - 1.0) ** 2) ** 0.5
        return _update_history(result, "KGE")
    else:
        if axis is None:
            from scipy.stats import pearsonr

            obsc, modc = matchedcompressed(obs, mod)
            if len(obsc) < 2:
                r = 0.0
            else:
                r, _ = pearsonr(obsc, modc)
        else:
            # Manual vectorized correlation for numpy with axis
            obs_mean = np.nanmean(obs, axis=axis, keepdims=True)
            mod_mean = np.nanmean(mod, axis=axis, keepdims=True)
            obs_std = obs - obs_mean
            mod_std = mod - mod_mean
            num = np.nansum(obs_std * mod_std, axis=axis)
            den = np.sqrt(np.nansum(obs_std**2, axis=axis) * np.nansum(mod_std**2, axis=axis))
            with np.errstate(divide="ignore", invalid="ignore"):
                r = num / den
                r = np.where(np.isnan(r), 0.0, r)

        alpha = np.ma.std(mod, axis=axis) / np.ma.std(obs, axis=axis)
        beta = np.ma.mean(mod, axis=axis) / np.ma.mean(obs, axis=axis)
        result = 1.0 - ((r - 1.0) ** 2 + (alpha - 1.0) ** 2 + (beta - 1.0) ** 2) ** 0.5
        return result.item() if np.ndim(result) == 0 else result

MAE(obs, mod, axis=None)

Mean Absolute Error (MAE).

Typical Use Cases

  • Quantifying the average magnitude of errors between model and observations, regardless of direction.
  • Used in model evaluation, forecast verification, and regression analysis.

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values. mod : numpy.ndarray or xarray.DataArray Model or predicted values. axis : int, str, or iterable of such, optional Axis or dimension along which to compute MAE.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Mean absolute error.

Examples

import numpy as np from monet_stats.error_metrics import MAE obs = np.array([1, 2, 3]) mod = np.array([2, 2, 4]) MAE(obs, mod) 0.6666666666666666

Source code in src/monet_stats/error_metrics.py
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
def MAE(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Mean Absolute Error (MAE).

    Typical Use Cases
    -----------------
    - Quantifying the average magnitude of errors between model and observations,
      regardless of direction.
    - Used in model evaluation, forecast verification, and regression analysis.

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values.
    mod : numpy.ndarray or xarray.DataArray
        Model or predicted values.
    axis : int, str, or iterable of such, optional
        Axis or dimension along which to compute MAE.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Mean absolute error.

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.error_metrics import MAE
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> MAE(obs, mod)
    0.6666666666666666
    """
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        dim = _resolve_axis_to_dim(obs, axis)
        result = abs(mod - obs).mean(dim=dim, keep_attrs=True)
        return _update_history(result, "MAE")
    else:
        result = np.ma.abs(np.subtract(mod, obs)).mean(axis=axis)
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

MAPE(obs, mod, axis=None)

Mean Absolute Percentage Error (MAPE).

Typical Use Cases

  • Quantifying the average relative error between model and observations as a percentage.
  • Used in time series forecasting, regression, and model evaluation for percentage-based error assessment.

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values. mod : numpy.ndarray or xarray.DataArray Model or predicted values. axis : int, str, or iterable of such, optional Axis or dimension along which to compute MAPE.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Mean absolute percentage error (in percent).

Examples

import numpy as np from monet_stats.error_metrics import MAPE obs = np.array([1, 2, 3]) mod = np.array([2, 2, 4]) MAPE(obs, mod) 50.0

Source code in src/monet_stats/error_metrics.py
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
def MAPE(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Mean Absolute Percentage Error (MAPE).

    Typical Use Cases
    -----------------
    - Quantifying the average relative error between model and observations
      as a percentage.
    - Used in time series forecasting, regression, and model evaluation for
      percentage-based error assessment.

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values.
    mod : numpy.ndarray or xarray.DataArray
        Model or predicted values.
    axis : int, str, or iterable of such, optional
        Axis or dimension along which to compute MAPE.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Mean absolute percentage error (in percent).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.error_metrics import MAPE
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> MAPE(obs, mod)
    50.0
    """
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        dim = _resolve_axis_to_dim(obs, axis)
        result = (100 * abs(mod - obs) / abs(obs)).mean(dim=dim, keep_attrs=True)
        return _update_history(result, "MAPE")
    else:
        result = (100 * np.ma.abs(np.subtract(mod, obs)) / np.ma.abs(obs)).mean(axis=axis)
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

MASE(obs, mod, axis=None)

Mean Absolute Scaled Error (MASE).

Typical Use Cases

  • Quantifying model error relative to the error of a simple baseline model (e.g., naive forecast).
  • Used in time series forecasting and model evaluation.
  • Provides scale-independent comparison across different datasets.

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values. mod : numpy.ndarray or xarray.DataArray Model or predicted values. axis : int, str, or iterable of such, optional Axis or dimension along which to compute MASE.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Mean absolute scaled error (unitless).

Examples

import numpy as np from monet_stats.error_metrics import MASE obs = np.array([1, 2, 3, 4]) mod = np.array([1.1, 2.1, 3.1, 4.1]) MASE(obs, mod) 0.1

Source code in src/monet_stats/error_metrics.py
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
def MASE(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Mean Absolute Scaled Error (MASE).

    Typical Use Cases
    -----------------
    - Quantifying model error relative to the error of a simple baseline model
      (e.g., naive forecast).
    - Used in time series forecasting and model evaluation.
    - Provides scale-independent comparison across different datasets.

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values.
    mod : numpy.ndarray or xarray.DataArray
        Model or predicted values.
    axis : int, str, or iterable of such, optional
        Axis or dimension along which to compute MASE.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Mean absolute scaled error (unitless).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.error_metrics import MASE
    >>> obs = np.array([1, 2, 3, 4])
    >>> mod = np.array([1.1, 2.1, 3.1, 4.1])
    >>> MASE(obs, mod)
    0.1
    """
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        dim = _resolve_axis_to_dim(obs, axis)

        # Calculate naive forecast error (using previous observation)
        if "time" in obs.dims:
            naive_error = abs(obs - obs.shift(time=1)).mean(dim=dim, skipna=True)
        else:
            # Fallback if time is not named 'time'
            naive_error = abs(obs - obs.shift({obs.dims[0]: 1})).mean(dim=dim, skipna=True)

        model_error = abs(mod - obs).mean(dim=dim, keep_attrs=True)
        result = model_error / naive_error
        return _update_history(result, "MASE")
    else:
        # Calculate naive forecast error (using previous observation)
        if axis is not None:
            naive_diff = np.diff(obs, axis=axis)
            naive_error = np.mean(np.abs(naive_diff), axis=axis)
        else:
            naive_diff = np.diff(obs)
            naive_error = np.mean(np.abs(naive_diff))
        model_error = np.mean(np.abs(np.subtract(mod, obs)), axis=axis)
        result = model_error / naive_error
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

MSE(obs, mod, axis=None)

Mean Squared Error (MSE).

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values. mod : numpy.ndarray or xarray.DataArray Model predicted values. axis : int, str, or iterable of such, optional Axis or dimension along which to compute the error.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Mean squared error.

Examples

import numpy as np from monet_stats.error_metrics import MSE obs = np.array([1, 2, 3]) mod = np.array([2, 2, 4]) MSE(obs, mod) 0.6666666666666666

Source code in src/monet_stats/error_metrics.py
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
def MSE(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Mean Squared Error (MSE).

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values.
    mod : numpy.ndarray or xarray.DataArray
        Model predicted values.
    axis : int, str, or iterable of such, optional
        Axis or dimension along which to compute the error.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Mean squared error.

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.error_metrics import MSE
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> MSE(obs, mod)
    0.6666666666666666
    """
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        dim = _resolve_axis_to_dim(obs, axis)
        result = ((mod - obs) ** 2).mean(dim=dim, keep_attrs=True)
        return _update_history(result, "MSE")
    else:
        result = np.ma.mean((np.subtract(mod, obs)) ** 2, axis=axis)
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

NSE(obs, mod, axis=None)

Nash-Sutcliffe Efficiency (NSE).

Typical Use Cases

  • Quantifying the predictive power of hydrological models relative to the mean of observations.
  • Used in hydrology, meteorology, and environmental model evaluation.

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values. mod : numpy.ndarray or xarray.DataArray Model predicted values. axis : int, str, or iterable of such, optional Axis or dimension along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Nash-Sutcliffe efficiency (unitless).

Examples

import numpy as np from monet_stats.efficiency_metrics import NSE obs = np.array([1, 2, 3, 4]) mod = np.array([1.1, 2.1, 2.9, 4.1]) NSE(obs, mod) 0.992

Source code in src/monet_stats/efficiency_metrics.py
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
def NSE(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Nash-Sutcliffe Efficiency (NSE).

    Typical Use Cases
    -----------------
    - Quantifying the predictive power of hydrological models relative to the
      mean of observations.
    - Used in hydrology, meteorology, and environmental model evaluation.

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values.
    mod : numpy.ndarray or xarray.DataArray
        Model predicted values.
    axis : int, str, or iterable of such, optional
        Axis or dimension along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Nash-Sutcliffe efficiency (unitless).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.efficiency_metrics import NSE
    >>> obs = np.array([1, 2, 3, 4])
    >>> mod = np.array([1.1, 2.1, 2.9, 4.1])
    >>> NSE(obs, mod)
    0.992
    """
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        dim = _resolve_axis_to_dim(obs, axis)

        obs_mean = obs.mean(dim=dim)
        numerator = ((obs - mod) ** 2).sum(dim=dim)
        denominator = ((obs - obs_mean) ** 2).sum(dim=dim)

        # Handle division by zero
        result = 1.0 - (numerator / denominator)
        result = xr.where((numerator == 0) & (denominator == 0), 1.0, result)
        result = xr.where((numerator != 0) & (denominator == 0), -np.inf, result)

        return _update_history(result, "NSE")
    else:
        obs_mean = np.nanmean(obs, axis=axis, keepdims=True)
        numerator = np.nansum((obs - mod) ** 2, axis=axis)
        denominator = np.nansum((obs - obs_mean) ** 2, axis=axis)

        with np.errstate(divide="ignore", invalid="ignore"):
            result = 1.0 - (numerator / denominator)
            result = np.where((numerator == 0) & (denominator == 0), 1.0, result)
            result = np.where((numerator != 0) & (denominator == 0), -np.inf, result)
        return result.item() if np.ndim(result) == 0 else result

NSElog(obs, mod, axis=None)

Log Nash-Sutcliffe Efficiency (NSElog).

Calculates NSE on logarithmic-transformed data to focus on lower values.

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values (positive values only). mod : numpy.ndarray or xarray.DataArray Model predicted values (positive values only). axis : int, str, or iterable of such, optional Axis or dimension along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Log Nash-Sutcliffe efficiency (unitless).

Examples

import numpy as np from monet_stats.efficiency_metrics import NSElog obs = np.array([1, 10, 100]) mod = np.array([1.1, 9.0, 110]) NSElog(obs, mod) 0.988

Source code in src/monet_stats/efficiency_metrics.py
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
def NSElog(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Log Nash-Sutcliffe Efficiency (NSElog).

    Calculates NSE on logarithmic-transformed data to focus on lower values.

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values (positive values only).
    mod : numpy.ndarray or xarray.DataArray
        Model predicted values (positive values only).
    axis : int, str, or iterable of such, optional
        Axis or dimension along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Log Nash-Sutcliffe efficiency (unitless).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.efficiency_metrics import NSElog
    >>> obs = np.array([1, 10, 100])
    >>> mod = np.array([1.1, 9.0, 110])
    >>> NSElog(obs, mod)
    0.988
    """
    epsilon = 1e-6
    # Avoid double history update by using .data or being careful
    # We apply log first, then call NSE. NSE will handle history if it's a DataArray.
    obs_log = np.log(obs + epsilon)
    mod_log = np.log(mod + epsilon)
    result = NSE(obs_log, mod_log, axis=axis)

    # If it's a DataArray, NSE already updated history to "NSE".
    # We might want to "fix" it to "NSElog".
    if isinstance(result, xr.DataArray):
        # Update history specifically for NSElog
        return _update_history(result, "NSElog")
    return result

NSEm(obs, mod, axis=None)

Nash-Sutcliffe Efficiency (NSE) - robust to masked arrays.

This function is a wrapper for NSE that explicitly handles masked data and NaNs.

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values. mod : numpy.ndarray or xarray.DataArray Model predicted values. axis : int, str, or iterable of such, optional Axis or dimension along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Nash-Sutcliffe efficiency (unitless).

Examples

import numpy as np from monet_stats.efficiency_metrics import NSEm obs = np.array([1, 2, np.nan, 4]) mod = np.array([1.1, 2.1, 3.0, 4.1]) NSEm(obs, mod) 0.995

Source code in src/monet_stats/efficiency_metrics.py
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
def NSEm(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Nash-Sutcliffe Efficiency (NSE) - robust to masked arrays.

    This function is a wrapper for NSE that explicitly handles masked data and NaNs.

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values.
    mod : numpy.ndarray or xarray.DataArray
        Model predicted values.
    axis : int, str, or iterable of such, optional
        Axis or dimension along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Nash-Sutcliffe efficiency (unitless).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.efficiency_metrics import NSEm
    >>> obs = np.array([1, 2, np.nan, 4])
    >>> mod = np.array([1.1, 2.1, 3.0, 4.1])
    >>> NSEm(obs, mod)
    0.995
    """
    # Standard NSE implementation already handles NaNs if using nan-aware functions
    res = NSE(obs, mod, axis=axis)
    if isinstance(res, xr.DataArray):
        return _update_history(res, "NSEm")
    return res

PC(obs, mod, axis=None, tolerance=0.1)

Percent of Correct (PC).

Calculates the percentage of model predictions that are within a specified tolerance of the observations.

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values. mod : numpy.ndarray or xarray.DataArray Model predicted values. axis : int, str, or iterable of such, optional Axis or dimension along which to compute the statistic. tolerance : float, optional Fraction of observed value used as tolerance (default 0.1).

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Percent of correct predictions (0-100%).

Examples

import numpy as np from monet_stats.efficiency_metrics import PC obs = np.array([1, 2, 3, 4]) mod = np.array([1.05, 2.5, 2.95, 4.05]) PC(obs, mod) 75.0

Source code in src/monet_stats/efficiency_metrics.py
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
def PC(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
    tolerance: float = 0.1,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Percent of Correct (PC).

    Calculates the percentage of model predictions that are within a specified
    tolerance of the observations.

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values.
    mod : numpy.ndarray or xarray.DataArray
        Model predicted values.
    axis : int, str, or iterable of such, optional
        Axis or dimension along which to compute the statistic.
    tolerance : float, optional
        Fraction of observed value used as tolerance (default 0.1).

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Percent of correct predictions (0-100%).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.efficiency_metrics import PC
    >>> obs = np.array([1, 2, 3, 4])
    >>> mod = np.array([1.05, 2.5, 2.95, 4.05])
    >>> PC(obs, mod)
    75.0
    """
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        dim = _resolve_axis_to_dim(obs, axis)

        tol = tolerance * abs(obs)
        correct = abs(obs - mod) <= tol
        result = (correct.sum(dim=dim) / correct.count(dim=dim)) * 100.0

        return _update_history(result, "PC")
    else:
        obs = np.asanyarray(obs)
        mod = np.asanyarray(mod)
        mask = np.isnan(obs) | np.isnan(mod)

        tol = tolerance * np.abs(obs)
        correct = np.abs(obs - mod) <= tol

        total = np.sum(~mask, axis=axis)
        correct_sum = np.sum(correct & ~mask, axis=axis)

        with np.errstate(divide="ignore", invalid="ignore"):
            result = (correct_sum / total) * 100.0
            result = np.where(total == 0, np.nan, result)

        return result.item() if np.ndim(result) == 0 else result

mNSE(obs, mod, axis=None)

Modified Nash-Sutcliffe Efficiency (mNSE).

Uses absolute differences instead of squared differences.

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values. mod : numpy.ndarray or xarray.DataArray Model predicted values. axis : int, str, or iterable of such, optional Axis or dimension along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Modified Nash-Sutcliffe efficiency (unitless).

Examples

import numpy as np from monet_stats.efficiency_metrics import mNSE obs = np.array([1, 2, 3, 4]) mod = np.array([1.1, 2.1, 2.9, 4.1]) mNSE(obs, mod) 0.92

Source code in src/monet_stats/efficiency_metrics.py
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
def mNSE(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Modified Nash-Sutcliffe Efficiency (mNSE).

    Uses absolute differences instead of squared differences.

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values.
    mod : numpy.ndarray or xarray.DataArray
        Model predicted values.
    axis : int, str, or iterable of such, optional
        Axis or dimension along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Modified Nash-Sutcliffe efficiency (unitless).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.efficiency_metrics import mNSE
    >>> obs = np.array([1, 2, 3, 4])
    >>> mod = np.array([1.1, 2.1, 2.9, 4.1])
    >>> mNSE(obs, mod)
    0.92
    """
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        dim = _resolve_axis_to_dim(obs, axis)

        obs_mean = obs.mean(dim=dim)
        numerator = abs(obs - mod).sum(dim=dim)
        denominator = abs(obs - obs_mean).sum(dim=dim)

        result = 1.0 - (numerator / denominator)
        result = xr.where((numerator == 0) & (denominator == 0), 1.0, result)
        result = xr.where((numerator != 0) & (denominator == 0), -np.inf, result)

        return _update_history(result, "mNSE")
    else:
        obs_mean = np.nanmean(obs, axis=axis, keepdims=True)
        numerator = np.nansum(np.abs(obs - mod), axis=axis)
        denominator = np.nansum(np.abs(obs - obs_mean), axis=axis)

        with np.errstate(divide="ignore", invalid="ignore"):
            result = 1.0 - (numerator / denominator)
            result = np.where((numerator == 0) & (denominator == 0), 1.0, result)
            result = np.where((numerator != 0) & (denominator == 0), -np.inf, result)
        return result.item() if np.ndim(result) == 0 else result

rNSE(obs, mod, axis=None)

Relative Nash-Sutcliffe Efficiency (rNSE).

Normalizes errors by the magnitude of observed values. Formula: 1 - [ sum( ((obs - mod)/obs)^2 ) / sum( ((obs - mean)/obs)^2 ) ]

Parameters

obs : numpy.ndarray or xarray.DataArray Observed values (should be non-zero for normalization). mod : numpy.ndarray or xarray.DataArray Model predicted values. axis : int, str, or iterable of such, optional Axis or dimension along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Relative Nash-Sutcliffe efficiency (unitless).

Examples

import numpy as np from monet_stats.efficiency_metrics import rNSE obs = np.array([1, 2, 3, 4]) mod = np.array([1.1, 2.1, 2.9, 4.1]) rNSE(obs, mod) 0.994261721483555

Source code in src/monet_stats/efficiency_metrics.py
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
def rNSE(
    obs: Union[np.ndarray, xr.DataArray],
    mod: Union[np.ndarray, xr.DataArray],
    axis: Optional[Union[int, str, Iterable[Union[int, str]]]] = None,
) -> Union[np.number, np.ndarray, xr.DataArray]:
    """
    Relative Nash-Sutcliffe Efficiency (rNSE).

    Normalizes errors by the magnitude of observed values.
    Formula: 1 - [ sum( ((obs - mod)/obs)^2 ) / sum( ((obs - mean)/obs)^2 ) ]

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed values (should be non-zero for normalization).
    mod : numpy.ndarray or xarray.DataArray
        Model predicted values.
    axis : int, str, or iterable of such, optional
        Axis or dimension along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Relative Nash-Sutcliffe efficiency (unitless).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.efficiency_metrics import rNSE
    >>> obs = np.array([1, 2, 3, 4])
    >>> mod = np.array([1.1, 2.1, 2.9, 4.1])
    >>> rNSE(obs, mod)
    0.994261721483555
    """
    epsilon = 1e-8
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        dim = _resolve_axis_to_dim(obs, axis)

        obs_mean = obs.mean(dim=dim)
        obs_safe = xr.where(abs(obs) < epsilon, epsilon, obs)

        numerator = (((obs - mod) / obs_safe) ** 2).sum(dim=dim)
        denominator = (((obs - obs_mean) / obs_safe) ** 2).sum(dim=dim)

        result = 1.0 - (numerator / denominator)
        result = xr.where((numerator == 0) & (denominator == 0), 1.0, result)
        result = xr.where((numerator != 0) & (denominator == 0), -np.inf, result)

        return _update_history(result, "rNSE")
    else:
        obs_mean = np.nanmean(obs, axis=axis, keepdims=True)
        obs_safe = np.where(np.abs(obs) < epsilon, epsilon, obs)

        with np.errstate(divide="ignore", invalid="ignore"):
            numerator = np.nansum(((obs - mod) / obs_safe) ** 2, axis=axis)
            denominator = np.nansum(((obs - obs_mean) / obs_safe) ** 2, axis=axis)
            result = 1.0 - (numerator / denominator)
            result = np.where((numerator == 0) & (denominator == 0), 1.0, result)
            result = np.where((numerator != 0) & (denominator == 0), -np.inf, result)
        return result.item() if np.ndim(result) == 0 else result