Skip to content

Correlation Metrics

Statistical correlation and skill score calculations for model evaluation.

Correlation and Agreement Metrics for Model Evaluation

AC(obs, mod, axis=None)

Anomaly Correlation (AC).

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 along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Anomaly correlation coefficient (unitless, -1 to 1).

Examples

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

Source code in src/monet_stats/correlation_metrics.py
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
def AC(
    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]:
    """
    Anomaly Correlation (AC).

    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 along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Anomaly correlation coefficient (unitless, -1 to 1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import AC
    >>> obs = np.array([1, 2, 3, 4])
    >>> mod = np.array([1.1, 2.1, 2.9, 4.1])
    >>> AC(obs, mod)
    0.9922778767136677
    """
    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

        obs_bar = obs.mean(dim=dim)
        mod_bar = mod.mean(dim=dim)
        obs_anom = obs - obs_bar
        mod_anom = mod - mod_bar
        p1 = (mod_anom * obs_anom).sum(dim=dim)
        p2 = ((mod_anom**2).sum(dim=dim) * (obs_anom**2).sum(dim=dim)) ** 0.5
        result = p1 / p2
        return _update_history(result, "AC")
    else:
        obs_bar = np.ma.mean(obs, axis=axis)
        mod_bar = np.ma.mean(mod, axis=axis)
        if axis is not None:
            # Need to keep dims for subtraction if axis is not None
            obs_bar_kd = np.ma.mean(obs, axis=axis, keepdims=True)
            mod_bar_kd = np.ma.mean(mod, axis=axis, keepdims=True)
        else:
            obs_bar_kd = obs_bar
            mod_bar_kd = mod_bar
        obs_anom = np.subtract(obs, obs_bar_kd)
        mod_anom = np.subtract(mod, mod_bar_kd)
        p1 = np.ma.sum(np.ma.multiply(mod_anom, obs_anom), axis=axis)
        p2 = np.ma.sqrt(np.ma.multiply(np.ma.sum(obs_anom**2, axis=axis), np.ma.sum(mod_anom**2, axis=axis)))
        result = p1 / p2
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

CCC(obs, mod, axis=None)

Concordance Correlation Coefficient (CCC).

Typical Use Cases

  • Quantifying the agreement between model and observations, accounting for precision and accuracy.
  • Used in model evaluation to assess how well model predictions agree with observations.
  • Measures how far the values deviate from the line of perfect concordance (slope=1, intercept=0).

Typical Values and Range

  • Range: -1 to 1
  • 1: Perfect agreement between model and observations
  • 0: No agreement
  • -1: Perfect negative agreement

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 the coefficient.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Concordance correlation coefficient (unitless, -1 to 1).

Examples

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

Source code in src/monet_stats/correlation_metrics.py
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
def CCC(
    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]:
    """
    Concordance Correlation Coefficient (CCC).

    Typical Use Cases
    -----------------
    - Quantifying the agreement between model and observations, accounting for
      precision and accuracy.
    - Used in model evaluation to assess how well model predictions agree with
      observations.
    - Measures how far the values deviate from the line of perfect concordance
      (slope=1, intercept=0).

    Typical Values and Range
    ------------------------
    - Range: -1 to 1
    - 1: Perfect agreement between model and observations
    - 0: No agreement
    - -1: Perfect negative agreement

    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 the coefficient.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Concordance correlation coefficient (unitless, -1 to 1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import CCC
    >>> obs = np.array([1, 2, 3, 4])
    >>> mod = np.array([1.1, 2.1, 2.9, 4.1])
    >>> CCC(obs, mod)
    0.9984779299847792
    """
    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

        # Calculate means
        obs_mean = obs.mean(dim=dim)
        mod_mean = mod.mean(dim=dim)

        # Calculate variances and covariance
        obs_var = obs.var(dim=dim)
        mod_var = mod.var(dim=dim)
        covar = ((obs - obs_mean) * (mod - mod_mean)).mean(dim=dim)

        # Calculate CCC
        numerator = 2 * covar
        denominator = obs_var + mod_var + (obs_mean - mod_mean) ** 2
        result = numerator / denominator
        return _update_history(result, "CCC")
    else:
        # Calculate means
        obs_mean = np.nanmean(obs, axis=axis)
        mod_mean = np.nanmean(mod, axis=axis)

        # Calculate variances and covariance
        obs_var = np.nanvar(obs, axis=axis)
        mod_var = np.nanvar(mod, axis=axis)
        if axis is not None:
            obs_mean_kd = np.nanmean(obs, axis=axis, keepdims=True)
            mod_mean_kd = np.nanmean(mod, axis=axis, keepdims=True)
        else:
            obs_mean_kd = obs_mean
            mod_mean_kd = mod_mean
        covar = np.nanmean((obs - obs_mean_kd) * (mod - mod_mean_kd), axis=axis)

        # Calculate CCC
        numerator = 2 * covar
        denominator = obs_var + mod_var + (obs_mean - mod_mean) ** 2
        result = numerator / denominator
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

E1(obs, mod, axis=None)

Modified Coefficient of Efficiency (E1).

Typical Use Cases

  • Quantifying the efficiency of model predictions relative to observed mean, robust to outliers.
  • Used in hydrology, meteorology, and model skill assessment.

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 along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Modified coefficient of efficiency (unitless, -inf to 1).

Examples

import numpy as np from monet_stats.correlation_metrics import E1 obs = np.array([1, 2, 3]) mod = np.array([2, 2, 4]) E1(obs, mod) 0.0

Source code in src/monet_stats/correlation_metrics.py
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
def E1(
    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 Coefficient of Efficiency (E1).

    Typical Use Cases
    -----------------
    - Quantifying the efficiency of model predictions relative to observed mean,
      robust to outliers.
    - Used in hydrology, meteorology, and model skill assessment.

    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 along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Modified coefficient of efficiency (unitless, -inf to 1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import E1
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> E1(obs, mod)
    0.0
    """
    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

        num = abs(obs - mod).sum(dim=dim)
        denom = abs(obs - obs.mean(dim=dim)).sum(dim=dim)
        result = 1.0 - (num / denom)
        result = xr.where((num == 0) & (denom == 0), 1.0, result)
        result = xr.where((num != 0) & (denom == 0), -np.inf, result)

        return _update_history(result, "E1")
    else:
        num = np.ma.abs(np.subtract(obs, mod)).sum(axis=axis)
        mean_obs = np.ma.mean(obs, axis=axis, keepdims=True)
        denom = np.ma.abs(np.subtract(obs, mean_obs)).sum(axis=axis)
        with np.errstate(divide="ignore", invalid="ignore"):
            result = 1.0 - (num / denom)
            result = np.where((num == 0) & (denom == 0), 1.0, result)
            result = np.where((num != 0) & (denom == 0), -np.inf, result)
        return result.item() if np.ndim(result) == 0 else result

E1_prime(obs, mod, axis=None)

Modified Coefficient of Efficiency (E1') - Alternative formulation.

Typical Use Cases

  • Quantifying the efficiency of model predictions relative to observed mean, robust to outliers.
  • Used in hydrology, meteorology, and model skill assessment as an alternative to E1.

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 along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Modified coefficient of efficiency (unitless, -inf to 1).

Examples

import numpy as np from monet_stats.correlation_metrics import E1_prime obs = np.array([1, 2, 3]) mod = np.array([2, 2, 4]) E1_prime(obs, mod) 0.0

Source code in src/monet_stats/correlation_metrics.py
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
def E1_prime(
    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 Coefficient of Efficiency (E1') - Alternative formulation.

    Typical Use Cases
    -----------------
    - Quantifying the efficiency of model predictions relative to observed mean,
      robust to outliers.
    - Used in hydrology, meteorology, and model skill assessment as an
      alternative to E1.

    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 along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Modified coefficient of efficiency (unitless, -inf to 1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import E1_prime
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> E1_prime(obs, mod)
    0.0
    """
    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

        obs_mean = obs.mean(dim=dim)
        num = abs(obs - mod).sum(dim=dim)
        denom = abs(obs - obs_mean).sum(dim=dim)
        # Handle case where denominator is 0
        result = 1.0 - (num / denom)
        result = xr.where((num == 0) & (denom == 0), 1.0, result)
        result = xr.where((num != 0) & (denom == 0), -np.inf, result)

        return _update_history(result, "E1_prime")
    else:
        if axis is None:
            obs_c, mod_c = matchedcompressed(obs, mod)
            obs_mean_kd = np.nanmean(obs_c)
        else:
            obs_c, mod_c = obs, mod
            obs_mean_kd = np.nanmean(obs_c, axis=axis, keepdims=True)

        num = np.nansum(np.abs(obs_c - mod_c), axis=axis)
        denom = np.nansum(np.abs(obs_c - obs_mean_kd), axis=axis)
        with np.errstate(divide="ignore", invalid="ignore"):
            result = 1.0 - (num / denom)
            if np.ndim(result) == 0:
                if num == 0 and denom == 0:
                    result = np.array(1.0)
                elif denom == 0:
                    result = np.array(-np.inf)
            else:
                result = np.where((num == 0) & (denom == 0), 1.0, result)
                result = np.where((num != 0) & (denom == 0), -np.inf, result)
        return result.item() if np.ndim(result) == 0 else result

IOA(obs, mod, axis=None)

Index of Agreement (IOA).

Typical Use Cases

  • Quantifying the agreement between model and observations, normalized by total deviation.
  • Used in model evaluation for skill 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 IOA.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Index of agreement (unitless, 0-1).

Examples

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

Source code in src/monet_stats/error_metrics.py
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
def IOA(
    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]:
    """
    Index of Agreement (IOA).

    Typical Use Cases
    -----------------
    - Quantifying the agreement between model and observations, normalized by
      total deviation.
    - Used in model evaluation for skill 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 IOA.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Index of agreement (unitless, 0-1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.error_metrics import IOA
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> IOA(obs, mod)
    0.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)
        num = ((obs - mod) ** 2).sum(dim=dim)
        denom = ((abs(mod - obs_mean) + abs(obs - obs_mean)) ** 2).sum(dim=dim)
        result = 1.0 - (num / denom)
        return _update_history(result, "IOA")
    else:
        obs_m = np.ma.masked_invalid(obs)
        mod_m = np.ma.masked_invalid(mod)
        obs_mean = np.ma.mean(obs_m, axis=axis, keepdims=True)
        num = np.ma.sum((obs_m - mod_m) ** 2, axis=axis)
        denom = np.ma.sum((np.ma.abs(mod_m - obs_mean) + np.ma.abs(obs_m - obs_mean)) ** 2, axis=axis)
        with np.errstate(divide="ignore", invalid="ignore"):
            result = 1.0 - (num / denom)
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

IOA_prime(obs, mod, axis=None)

Index of Agreement (IOA') - Alternative formulation.

Typical Use Cases

  • Quantifying the agreement between model and observations, normalized by total deviation.
  • Used in model evaluation for skill assessment as an alternative to IOA.

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 along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Index of agreement (unitless, 0-1).

Examples

import numpy as np from monet_stats.correlation_metrics import IOA_prime obs = np.array([1, 2, 3]) mod = np.array([2, 2, 4]) IOA_prime(obs, mod) 0.8

Source code in src/monet_stats/correlation_metrics.py
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
def IOA_prime(
    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]:
    """
    Index of Agreement (IOA') - Alternative formulation.

    Typical Use Cases
    -----------------
    - Quantifying the agreement between model and observations, normalized by
      total deviation.
    - Used in model evaluation for skill assessment as an alternative to IOA.

    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 along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Index of agreement (unitless, 0-1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import IOA_prime
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> IOA_prime(obs, mod)
    0.8
    """
    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

        obsmean = obs.mean(dim=dim)
        num = ((obs - mod) ** 2).sum(dim=dim)
        denom = ((abs(mod - obsmean) + abs(obs - obsmean)) ** 2).sum(dim=dim)
        # Handle case where denominator is 0
        result = 1.0 - (num / denom)
        result = xr.where((num == 0) & (denom == 0), 1.0, result)
        result = xr.where((num != 0) & (denom == 0), -np.inf, result)

        return _update_history(result, "IOA_prime")
    else:
        if axis is None:
            obs_c, mod_c = matchedcompressed(obs, mod)
            obsmean_kd = np.nanmean(obs_c)
        else:
            obs_c, mod_c = obs, mod
            obsmean_kd = np.nanmean(obs_c, axis=axis, keepdims=True)

        num = np.nansum((obs_c - mod_c) ** 2, axis=axis)
        denom = np.nansum((np.abs(mod_c - obsmean_kd) + np.abs(obs_c - obsmean_kd)) ** 2, axis=axis)
        with np.errstate(divide="ignore", invalid="ignore"):
            result = 1.0 - (num / denom)
            if np.ndim(result) == 0:
                if num == 0 and denom == 0:
                    result = np.array(1.0)
                elif denom == 0:
                    result = np.array(-np.inf)
            else:
                result = np.where((num == 0) & (denom == 0), 1.0, result)
                result = np.where((num != 0) & (denom == 0), -np.inf, result)
        return result.item() if np.ndim(result) == 0 else result

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

R2(obs, mod, axis=None)

Coefficient of Determination (R^2, unitless).

Typical Use Cases

  • Quantifying how well model predictions explain the variance in observations.
  • Used in regression analysis, model skill assessment, and forecast verification.

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 Coefficient of determination (R^2).

Examples

import numpy as np from monet_stats.correlation_metrics import R2 obs = np.array([1, 2, 3, 4]) mod = np.array([1.1, 1.9, 3.2, 3.8]) R2(obs, mod) 0.9846153846153847

Source code in src/monet_stats/correlation_metrics.py
 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
 78
 79
 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
def R2(
    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]:
    """
    Coefficient of Determination (R^2, unitless).

    Typical Use Cases
    -----------------
    - Quantifying how well model predictions explain the variance in observations.
    - Used in regression analysis, model skill assessment, and forecast
      verification.

    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
        Coefficient of determination (R^2).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import R2
    >>> obs = np.array([1, 2, 3, 4])
    >>> mod = np.array([1.1, 1.9, 3.2, 3.8])
    >>> R2(obs, mod)
    0.9846153846153847
    """
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        if axis is None:
            # Default to all dimensions if None
            dim = obs.dims
        elif isinstance(axis, int):
            dim = obs.dims[axis]
        else:
            dim = axis

        # Use native xarray correlation for speed and laziness (Aero Protocol)
        r = xr.corr(obs, mod, dim=dim)
        result = r**2
        return _update_history(result, "R2")
    else:
        from scipy.stats import pearsonr

        if axis is None:
            obsc, modc = matchedcompressed(obs, mod)
            if len(obsc) < 2 or np.var(obsc) == 0 or np.var(modc) == 0:
                return 0.0
            r_val, _ = pearsonr(obsc, modc)
            if np.isnan(r_val):
                return 0.0
            return r_val**2
        else:
            # Manual vectorized R2
            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
                result = np.where(np.isnan(r), 0.0, r**2)
                return result.item() if np.ndim(result) == 0 else result

RMSE(obs, mod, axis=None)

Root Mean Square Error (RMSE).

Typical Use Cases

  • Quantifying the average magnitude of errors between model and observations, accounting for large errors more heavily than MAE.
  • 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 RMSE.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Root mean square error.

Examples

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

Source code in src/monet_stats/error_metrics.py
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
def RMSE(
    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]:
    """
    Root Mean Square Error (RMSE).

    Typical Use Cases
    -----------------
    - Quantifying the average magnitude of errors between model and observations,
      accounting for large errors more heavily than MAE.
    - 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 RMSE.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Root mean square error.

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.error_metrics import RMSE
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> RMSE(obs, mod)
    0.816496580927726
    """
    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) ** 0.5
        return _update_history(result, "RMSE")
    else:
        result = np.ma.sqrt(np.ma.mean((np.subtract(mod, obs)) ** 2, axis=axis))
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

RMSEs(obs, mod, axis=None)

Root Mean Squared Error between observations and regression fit.

(RMSEs, model unit)

Typical Use Cases

  • Quantifying the error between observations and a regression fit to the model predictions.
  • Used in model evaluation to assess how well a regression fit to the model matches 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 along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray, optional Root mean squared error value(s).

Examples

import numpy as np from monet_stats.correlation_metrics import RMSEs obs = np.array([1, 2, 3, 4]) mod = np.array([2, 2, 2, 2]) RMSEs(obs, mod) 0.7071067811865476

Source code in src/monet_stats/correlation_metrics.py
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
def RMSEs(
    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, None]:
    """
    Root Mean Squared Error between observations and regression fit.

    (RMSEs, model unit)

    Typical Use Cases
    -----------------
    - Quantifying the error between observations and a regression fit to the
      model predictions.
    - Used in model evaluation to assess how well a regression fit to the model
      matches 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 along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray, optional
        Root mean squared error value(s).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import RMSEs
    >>> obs = np.array([1, 2, 3, 4])
    >>> mod = np.array([2, 2, 2, 2])
    >>> RMSEs(obs, mod)
    0.7071067811865476
    """
    res = _vectorized_regression_stats(obs, mod, axis=axis, mode="RMSEs")
    return _update_history(res, "RMSEs")

RMSEu(obs, mod, axis=None)

Root Mean Squared Error between regression fit and model predictions.

(RMSEu, model unit)

Typical Use Cases

  • Quantifying the error between a linear regression fit to observations and the model predictions.
  • Used in model evaluation to assess how well a regression fit to obs matches the model output.

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 along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray, optional Root mean squared error value(s).

Examples

import numpy as np from monet_stats.correlation_metrics import RMSEu obs = np.array([1, 2, 3, 4]) mod = np.array([2, 2, 2, 2]) RMSEu(obs, mod) 0.7071067811865476

Source code in src/monet_stats/correlation_metrics.py
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
def RMSEu(
    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, None]:
    """
    Root Mean Squared Error between regression fit and model predictions.

    (RMSEu, model unit)

    Typical Use Cases
    -----------------
    - Quantifying the error between a linear regression fit to observations and
      the model predictions.
    - Used in model evaluation to assess how well a regression fit to obs
      matches the model output.

    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 along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray, optional
        Root mean squared error value(s).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import RMSEu
    >>> obs = np.array([1, 2, 3, 4])
    >>> mod = np.array([2, 2, 2, 2])
    >>> RMSEu(obs, mod)
    0.7071067811865476
    """
    res = _vectorized_regression_stats(obs, mod, axis=axis, mode="RMSEu")
    return _update_history(res, "RMSEu")

WDAC(obs, mod, axis=None)

Wind Direction Anomaly Correlation (WDAC).

Parameters

obs : numpy.ndarray or xarray.DataArray Observed wind direction values (degrees). mod : numpy.ndarray or xarray.DataArray Modeled wind direction values (degrees). axis : int, str, or iterable of such, optional Axis along which to compute the metric.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray WDAC value(s).

Examples

import numpy as np from monet_stats.correlation_metrics import WDAC obs = np.array([350, 10, 20]) mod = np.array([345, 15, 25]) WDAC(obs, mod) 0.9992386127814763

Source code in src/monet_stats/correlation_metrics.py
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
def WDAC(
    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]:
    """
    Wind Direction Anomaly Correlation (WDAC).

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed wind direction values (degrees).
    mod : numpy.ndarray or xarray.DataArray
        Modeled wind direction values (degrees).
    axis : int, str, or iterable of such, optional
        Axis along which to compute the metric.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        WDAC value(s).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import WDAC
    >>> obs = np.array([350, 10, 20])
    >>> mod = np.array([345, 15, 25])
    >>> WDAC(obs, mod)
    0.9992386127814763
    """
    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

        obs_rad = obs * np.pi / 180.0
        mod_rad = mod * np.pi / 180.0
        obs_anom = obs_rad - obs_rad.mean(dim=dim)
        mod_anom = mod_rad - mod_rad.mean(dim=dim)
        numerator = (np.sin(obs_anom) * np.sin(mod_anom)).sum(dim=dim)
        denominator = np.sqrt((np.sin(obs_anom) ** 2).sum(dim=dim) * (np.sin(mod_anom) ** 2).sum(dim=dim))
        result = numerator / denominator
        return _update_history(result, "WDAC")
    else:
        obs_rad = np.deg2rad(obs)
        mod_rad = np.deg2rad(mod)
        if axis is not None:
            obs_bar_rad = np.ma.mean(obs_rad, axis=axis, keepdims=True)
            mod_bar_rad = np.ma.mean(mod_rad, axis=axis, keepdims=True)
        else:
            obs_bar_rad = np.ma.mean(obs_rad)
            mod_bar_rad = np.ma.mean(mod_rad)

        obs_anom = obs_rad - obs_bar_rad
        mod_anom = mod_rad - mod_bar_rad
        numerator = np.ma.sum(np.sin(obs_anom) * np.sin(mod_anom), axis=axis)
        denominator = np.ma.sqrt(
            np.ma.sum(np.sin(obs_anom) ** 2, axis=axis) * np.ma.sum(np.sin(mod_anom) ** 2, axis=axis)
        )
        result = numerator / denominator
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

WDIOA(obs, mod, axis=None)

Wind Direction Index of Agreement (WDIOA).

Standard version.

Typical Use Cases

  • Quantifying the agreement between observed and modeled wind directions, accounting for circularity.
  • Used in wind energy, meteorology, and air quality studies to assess wind direction model performance.

Typical Values and Range

  • Range: 0 to 1
  • 1: Perfect agreement between observed and modeled wind directions
  • 0: No agreement (as bad as using the mean of observations)

Parameters

obs : numpy.ndarray or xarray.DataArray Observed wind direction values (degrees). mod : numpy.ndarray or xarray.DataArray Modeled wind direction values (degrees). axis : int, str, or iterable of such, optional Axis along which to compute the metric.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Wind direction index of agreement (unitless, 0-1).

Examples

import numpy as np from monet_stats.correlation_metrics import WDIOA obs = np.array([350, 10, 20]) mod = np.array([345, 15, 25]) WDIOA(obs, mod) 0.8

Source code in src/monet_stats/correlation_metrics.py
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
def WDIOA(
    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]:
    """
    Wind Direction Index of Agreement (WDIOA).

    Standard version.

    Typical Use Cases
    -----------------
    - Quantifying the agreement between observed and modeled wind directions,
      accounting for circularity.
    - Used in wind energy, meteorology, and air quality studies to assess wind
      direction model performance.

    Typical Values and Range
    ------------------------
    - Range: 0 to 1
    - 1: Perfect agreement between observed and modeled wind directions
    - 0: No agreement (as bad as using the mean of observations)

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed wind direction values (degrees).
    mod : numpy.ndarray or xarray.DataArray
        Modeled wind direction values (degrees).
    axis : int, str, or iterable of such, optional
        Axis along which to compute the metric.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Wind direction index of agreement (unitless, 0-1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import WDIOA
    >>> obs = np.array([350, 10, 20])
    >>> mod = np.array([345, 15, 25])
    >>> WDIOA(obs, mod)
    0.8
    """
    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

        num = abs(circlebias(obs - mod)).sum(dim=dim)
        mean_obs = obs.mean(dim=dim)
        denom = (abs(circlebias(mod - mean_obs)) + abs(circlebias(obs - mean_obs))).sum(dim=dim)

        result = 1.0 - (num / denom)
        result = xr.where(denom == 0, 1.0, result)

        return _update_history(result, "WDIOA")
    else:
        num = np.ma.sum(np.ma.abs(circlebias(np.subtract(obs, mod))), axis=axis)
        mean_obs = np.ma.mean(obs, axis=axis, keepdims=True)
        denom = np.ma.sum(
            np.ma.abs(circlebias(np.subtract(mod, mean_obs))) + np.ma.abs(circlebias(np.subtract(obs, mean_obs))),
            axis=axis,
        )
        result = np.where(denom == 0, 1.0, 1.0 - (num / denom))
        return result.item() if np.ndim(result) == 0 else result

WDIOA_m(obs, mod, axis=None)

Wind Direction Index of Agreement (WDIOA_m).

Robust to masked arrays.

Typical Use Cases

  • Quantifying the agreement between observed and modeled wind directions, accounting for circularity.
  • Used in wind energy, meteorology, and air quality studies to assess wind direction model performance.

Typical Values and Range

  • Range: 0 to 1
  • 1: Perfect agreement between observed and modeled wind directions
  • 0: No agreement (as bad as using the mean of observations)

Parameters

obs : numpy.ndarray or xarray.DataArray Observed wind direction values (degrees). mod : numpy.ndarray or xarray.DataArray Modeled wind direction values (degrees). axis : int, str, or iterable of such, optional Axis along which to compute the metric.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Wind direction index of agreement (unitless, 0-1).

Examples

import numpy as np from monet_stats.correlation_metrics import WDIOA_m obs = np.array([350, 10, 20]) mod = np.array([345, 15, 25]) WDIOA_m(obs, mod) 0.8

Source code in src/monet_stats/correlation_metrics.py
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
def WDIOA_m(
    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]:
    """
    Wind Direction Index of Agreement (WDIOA_m).

    Robust to masked arrays.

    Typical Use Cases
    -----------------
    - Quantifying the agreement between observed and modeled wind directions,
      accounting for circularity.
    - Used in wind energy, meteorology, and air quality studies to assess wind
      direction model performance.

    Typical Values and Range
    ------------------------
    - Range: 0 to 1
    - 1: Perfect agreement between observed and modeled wind directions
    - 0: No agreement (as bad as using the mean of observations)

    Parameters
    ----------
    obs : numpy.ndarray or xarray.DataArray
        Observed wind direction values (degrees).
    mod : numpy.ndarray or xarray.DataArray
        Modeled wind direction values (degrees).
    axis : int, str, or iterable of such, optional
        Axis along which to compute the metric.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Wind direction index of agreement (unitless, 0-1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import WDIOA_m
    >>> obs = np.array([350, 10, 20])
    >>> mod = np.array([345, 15, 25])
    >>> WDIOA_m(obs, mod)
    0.8
    """
    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

        obsmean = obs.mean(dim=dim)
        num = (abs(circlebias_m(obs - mod))).sum(dim=dim)
        denom = (abs(circlebias_m(mod - obsmean)) + abs(circlebias_m(obs - obsmean))).sum(dim=dim)

        result = 1.0 - (num / denom)
        result = xr.where(denom == 0, 1.0, result)

        return _update_history(result, "WDIOA_m")
    else:
        obsmean = np.ma.mean(obs, axis=axis, keepdims=True)
        num = np.ma.sum(np.ma.abs(circlebias_m(np.subtract(obs, mod))), axis=axis)
        denom = np.ma.sum(
            np.ma.abs(circlebias_m(np.subtract(mod, obsmean))) + np.ma.abs(circlebias_m(np.subtract(obs, obsmean))),
            axis=axis,
        )
        result = np.where(denom == 0, 1.0, 1.0 - (num / denom))
        return result.item() if np.ndim(result) == 0 else result

WDRMSE(obs, mod, axis=None)

Wind Direction Root Mean Square Error (WDRMSE, model unit).

Standard version.

Typical Use Cases

  • Quantifying the average magnitude of wind direction errors, accounting for circularity.
  • Used in wind energy, meteorology, and air quality studies to assess wind direction model performance.

Parameters

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

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Wind direction root mean square error (degrees).

Examples

import numpy as np from monet_stats.correlation_metrics import WDRMSE obs = np.array([350, 10, 20]) mod = np.array([10, 20, 30]) WDRMSE(obs, mod) 20.0

Source code in src/monet_stats/correlation_metrics.py
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
229
230
231
232
def WDRMSE(
    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]:
    """
    Wind Direction Root Mean Square Error (WDRMSE, model unit).

    Standard version.

    Typical Use Cases
    -----------------
    - Quantifying the average magnitude of wind direction errors, accounting for
      circularity.
    - Used in wind energy, meteorology, and air quality studies to assess wind
      direction model performance.

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

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Wind direction root mean square error (degrees).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import WDRMSE
    >>> obs = np.array([350, 10, 20])
    >>> mod = np.array([10, 20, 30])
    >>> WDRMSE(obs, mod)
    20.0
    """
    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

        result = (circlebias(mod - obs) ** 2).mean(dim=dim, keep_attrs=True) ** 0.5
        return _update_history(result, "WDRMSE")
    else:
        result = np.ma.sqrt(np.ma.mean((circlebias(np.subtract(mod, obs))) ** 2, axis=axis))
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

WDRMSE_m(obs, mod, axis=None)

Wind Direction Root Mean Square Error (WDRMSE, model unit).

Robust to masked arrays.

Typical Use Cases

  • Quantifying the average magnitude of wind direction errors, accounting for circularity, robust to masked arrays.
  • Used in wind energy, meteorology, and air quality studies to assess wind direction model performance.

Parameters

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

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Wind direction root mean square error (degrees).

Examples

import numpy as np from monet_stats.correlation_metrics import WDRMSE_m obs = np.array([350, 10, 20]) mod = np.array([10, 20, 30]) WDRMSE_m(obs, mod) 20.0

Source code in src/monet_stats/correlation_metrics.py
115
116
117
118
119
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
166
167
168
169
170
171
172
def WDRMSE_m(
    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]:
    """
    Wind Direction Root Mean Square Error (WDRMSE, model unit).

    Robust to masked arrays.

    Typical Use Cases
    -----------------
    - Quantifying the average magnitude of wind direction errors, accounting for
      circularity, robust to masked arrays.
    - Used in wind energy, meteorology, and air quality studies to assess wind
      direction model performance.

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

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Wind direction root mean square error (degrees).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import WDRMSE_m
    >>> obs = np.array([350, 10, 20])
    >>> mod = np.array([10, 20, 30])
    >>> WDRMSE_m(obs, mod)
    20.0
    """
    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

        result = (circlebias_m(mod - obs) ** 2).mean(dim=dim, keep_attrs=True) ** 0.5
        return _update_history(result, "WDRMSE_m")
    else:
        result = np.ma.sqrt(np.ma.mean((circlebias_m(np.subtract(mod, obs))) ** 2, axis=axis))
        return result.item() if hasattr(result, "item") and np.ndim(result) == 0 else result

d1(obs, mod, axis=None)

Modified Index of Agreement (d1).

Typical Use Cases

  • Quantifying the agreement between model and observations, less sensitive to outliers than IOA.
  • Used in model evaluation for robust skill assessment.

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 along which to compute the statistic.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Modified index of agreement (unitless, 0-1).

Examples

import numpy as np from monet_stats.correlation_metrics import d1 obs = np.array([1, 2, 3]) mod = np.array([2, 2, 4]) d1(obs, mod) 0.5

Source code in src/monet_stats/correlation_metrics.py
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
def d1(
    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 Index of Agreement (d1).

    Typical Use Cases
    -----------------
    - Quantifying the agreement between model and observations, less sensitive
      to outliers than IOA.
    - Used in model evaluation for robust skill assessment.

    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 along which to compute the statistic.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Modified index of agreement (unitless, 0-1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import d1
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> d1(obs, mod)
    0.5
    """
    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

        num = abs(obs - mod).sum(dim=dim)
        mean_obs = obs.mean(dim=dim)
        denom = (abs(mod - mean_obs) + abs(obs - mean_obs)).sum(dim=dim)
        result = 1.0 - (num / denom)
        result = xr.where((num == 0) & (denom == 0), 1.0, result)
        result = xr.where((num != 0) & (denom == 0), -np.inf, result)

        return _update_history(result, "d1")
    else:
        num = np.ma.abs(np.subtract(obs, mod)).sum(axis=axis)
        mean_obs = np.ma.mean(obs, axis=axis, keepdims=True)
        denom = (np.ma.abs(np.subtract(mod, mean_obs)) + np.ma.abs(np.subtract(obs, mean_obs))).sum(axis=axis)
        with np.errstate(divide="ignore", invalid="ignore"):
            result = 1.0 - (num / denom)
            result = np.where((num == 0) & (denom == 0), 1.0, result)
            result = np.where((num != 0) & (denom == 0), -np.inf, result)
        return result.item() if np.ndim(result) == 0 else result

kendalltau(obs, mod, axis=None)

Kendall tau correlation coefficient (Aero Protocol: Standardized).

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 name along which to compute the coefficient.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Kendall rank correlation coefficient.

Examples

import numpy as np from monet_stats.correlation_metrics import kendalltau obs = np.array([1, 2, 3]) mod = np.array([2, 2, 4]) kendalltau(obs, mod) 1.0

Source code in src/monet_stats/correlation_metrics.py
1178
1179
1180
1181
1182
1183
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
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
def kendalltau(
    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]:
    """
    Kendall tau correlation coefficient (Aero Protocol: Standardized).

    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 name along which to compute the coefficient.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Kendall rank correlation coefficient.

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import kendalltau
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> kendalltau(obs, mod)
    1.0
    """
    from scipy.stats import kendalltau as _kendalltau

    # Standardize to DataArray to leverage Xarray's apply_ufunc vectorization
    is_numpy = not isinstance(obs, xr.DataArray)
    if is_numpy:
        obs = xr.DataArray(obs)
        mod = xr.DataArray(mod)
        if axis is not None:
            if isinstance(axis, int):
                dim = obs.dims[axis]
            elif isinstance(axis, Iterable) and not isinstance(axis, str):
                dim = [obs.dims[i] for i in axis]
            else:
                dim = axis
        else:
            dim = obs.dims
    else:
        if isinstance(mod, xr.DataArray):
            obs, mod = xr.align(obs, mod, join="inner")
        if axis is None:
            dim = obs.dims
        elif isinstance(axis, int):
            dim = obs.dims[axis]
        else:
            dim = axis

    if axis is None:
        obsc, modc = matchedcompressed(obs, mod)
        if len(obsc) < 2:
            return np.nan
        return _kendalltau(obsc, modc)[0]

    def _kendalltau_onlytau(a, b):
        a_flat = a.ravel()
        b_flat = b.ravel()
        mask = ~np.isnan(a_flat) & ~np.isnan(b_flat)
        if np.sum(mask) < 2:
            return np.nan
        return _kendalltau(a_flat[mask], b_flat[mask])[0]

    # Use apply_ufunc to eliminate manual loops and support Dask
    icd = [dim] if isinstance(dim, (str, int)) else list(dim)
    result = xr.apply_ufunc(
        _kendalltau_onlytau,
        obs,
        mod,
        input_core_dims=[icd] * 2,
        output_core_dims=[[]],
        vectorize=True,
        dask="parallelized",
        dask_gufunc_kwargs={"allow_rechunk": True},
        output_dtypes=[float],
    )

    if is_numpy:
        return result.values
    return _update_history(result, "kendalltau")

pearsonr(obs, mod, axis=None)

Pearson correlation coefficient.

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 name along which to compute the coefficient.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Pearson correlation coefficient.

Examples

import numpy as np from monet_stats.correlation_metrics import pearsonr obs = np.array([1, 2, 3]) mod = np.array([2, 4, 6]) pearsonr(obs, mod) 1.0

Source code in src/monet_stats/correlation_metrics.py
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
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
def pearsonr(
    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]:
    """
    Pearson correlation coefficient.

    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 name along which to compute the coefficient.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Pearson correlation coefficient.

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import pearsonr
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 4, 6])
    >>> pearsonr(obs, mod)
    1.0
    """
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        if axis is None:
            dim = obs.dims
        elif isinstance(axis, int):
            dim = obs.dims[axis]
        else:
            dim = axis

        # Use native xarray correlation for speed and laziness (Aero Protocol)
        result = xr.corr(obs, mod, dim=dim)
        return _update_history(result, "pearsonr")
    else:
        from scipy.stats import pearsonr as _pearsonr

        if axis is None:
            obsc, modc = matchedcompressed(obs, mod)
            if len(obsc) < 2 or np.var(obsc) == 0 or np.var(modc) == 0:
                return 0.0
            r_val, _ = _pearsonr(obsc, modc)
            return r_val if not np.isnan(r_val) else 0.0
        else:
            # For numpy with axis, use manual vectorized correlation with pairwise deletion
            obs = np.asanyarray(obs)
            mod = np.asanyarray(mod)
            mask = np.isnan(obs) | np.isnan(mod)
            obs = np.where(mask, np.nan, obs)
            mod = np.where(mask, np.nan, mod)

            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"):
                result = num / den
                return result.item() if np.ndim(result) == 0 else result

spearmanr(obs, mod, axis=None)

Spearman rank correlation coefficient (Aero Protocol: Vectorized).

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 the coefficient.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Spearman rank correlation coefficient.

Examples

import numpy as np from monet_stats.correlation_metrics import spearmanr obs = np.array([1, 2, 3]) mod = np.array([2, 2, 4]) spearmanr(obs, mod) 0.8660254037844387

Source code in src/monet_stats/correlation_metrics.py
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
def spearmanr(
    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]:
    """
    Spearman rank correlation coefficient (Aero Protocol: Vectorized).

    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 the coefficient.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Spearman rank correlation coefficient.

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import spearmanr
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([2, 2, 4])
    >>> spearmanr(obs, mod)
    0.8660254037844387
    """
    # Handle Xarray/NumPy alignment and dimension resolution
    if isinstance(obs, xr.DataArray) and isinstance(mod, xr.DataArray):
        obs, mod = xr.align(obs, mod, join="inner")
        if axis is None:
            dim = obs.dims
        elif isinstance(axis, int):
            dim = obs.dims[axis]
        else:
            dim = axis
    else:
        dim = axis

    if axis is None:
        obsc, modc = matchedcompressed(obs, mod)
        if len(obsc) < 2:
            return np.nan
        from scipy.stats import spearmanr as _spearmanr

        return _spearmanr(obsc, modc)[0]

    # Spearman is Pearson on ranks. Use rankdata along the axis.
    from scipy.stats import rankdata

    # Apply pairwise masking to ensure ranks are computed on the same set of points
    mask = np.isnan(obs) | np.isnan(mod)
    obs = xr.where(mask, np.nan, obs) if isinstance(obs, xr.DataArray) else np.where(mask, np.nan, obs)
    mod = xr.where(mask, np.nan, mod) if isinstance(mod, xr.DataArray) else np.where(mask, np.nan, mod)

    def _rank_wrapper(data, axis):
        # Handle all-NaN slices to avoid Scipy warnings or errors
        if np.all(np.isnan(data)):
            return np.full(data.shape, np.nan)
        return rankdata(data, axis=axis, nan_policy="omit")

    # Use apply_ufunc for rankdata to support both NumPy and Xarray/Dask
    if isinstance(obs, xr.DataArray):
        icd = [dim] if isinstance(dim, (str, int)) else list(dim)
        obs_ranked = xr.apply_ufunc(
            _rank_wrapper,
            obs,
            input_core_dims=[icd],
            output_core_dims=[icd],
            kwargs={"axis": -1},
            dask="parallelized",
            dask_gufunc_kwargs={"allow_rechunk": True},
            output_dtypes=[float],
            keep_attrs=True,
        )
        mod_ranked = xr.apply_ufunc(
            _rank_wrapper,
            mod,
            input_core_dims=[icd],
            output_core_dims=[icd],
            kwargs={"axis": -1},
            dask="parallelized",
            dask_gufunc_kwargs={"allow_rechunk": True},
            output_dtypes=[float],
            keep_attrs=True,
        )
        return pearsonr(obs_ranked, mod_ranked, axis=dim)
    else:
        obs_ranked = _rank_wrapper(obs, axis=axis)
        mod_ranked = _rank_wrapper(mod, axis=axis)
        return pearsonr(obs_ranked, mod_ranked, axis=axis)

taylor_skill(obs, mod, axis=None)

Taylor Skill Score (TSS).

Typical Use Cases

  • Summarizing model performance in a single skill score for use in Taylor diagrams.
  • Used in climate, weather, and environmental model evaluation.

Typical Values and Range

  • Range: 0 to 1
  • 1: Perfect agreement between model and observations
  • 0: No 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 the skill score.

Returns

numpy.number, numpy.ndarray, or xarray.DataArray Taylor skill score (unitless, 0-1).

Examples

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

Source code in src/monet_stats/correlation_metrics.py
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
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
881
882
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
def taylor_skill(
    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]:
    """
    Taylor Skill Score (TSS).

    Typical Use Cases
    -----------------
    - Summarizing model performance in a single skill score for use in Taylor
      diagrams.
    - Used in climate, weather, and environmental model evaluation.

    Typical Values and Range
    ------------------------
    - Range: 0 to 1
    - 1: Perfect agreement between model and observations
    - 0: No 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 the skill score.

    Returns
    -------
    numpy.number, numpy.ndarray, or xarray.DataArray
        Taylor skill score (unitless, 0-1).

    Examples
    --------
    >>> import numpy as np
    >>> from monet_stats.correlation_metrics import taylor_skill
    >>> obs = np.array([1, 2, 3])
    >>> mod = np.array([1.1, 1.9, 3.2])
    >>> taylor_skill(obs, mod)
    0.9995574044955781
    """
    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

        std_obs = obs.std(dim=dim)
        std_mod = mod.std(dim=dim)
        corr = xr.corr(obs, mod, dim=dim)

        # Calculate Taylor Skill Score using the common formula
        # S = 4 * (1 + R) / ( (sigma_p/sigma_o + sigma_o/sigma_p)^2 * (1 + R_max) )
        # Assuming R_max = 1.0
        norm_std = std_mod / std_obs
        result = (4.0 * (corr + 1.0)) / ((norm_std + 1.0 / norm_std) ** 2 * 2.0)
        return _update_history(result, "taylor_skill")
    else:
        std_obs = np.ma.std(obs, axis=axis)
        std_mod = np.ma.std(mod, axis=axis)
        from scipy.stats import pearsonr

        if axis is None:
            if np.ma.is_masked(obs):
                corr = pearsonr(obs.compressed(), mod.compressed())[0]
            else:
                corr = pearsonr(obs, mod)[0]
        else:
            # Vectorized correlation over axis for numpy
            obs_mean = np.nanmean(obs, axis=axis, keepdims=True)
            mod_mean = np.nanmean(mod, axis=axis, keepdims=True)
            obs_anom = obs - obs_mean
            mod_anom = mod - mod_mean
            num_corr = np.nansum(obs_anom * mod_anom, axis=axis)
            den_corr = np.sqrt(np.nansum(obs_anom**2, axis=axis) * np.nansum(mod_anom**2, axis=axis))
            with np.errstate(divide="ignore", invalid="ignore"):
                corr = num_corr / den_corr

        norm_std = std_mod / std_obs
        with np.errstate(divide="ignore", invalid="ignore"):
            result = (4.0 * (corr + 1.0)) / ((norm_std + 1.0 / norm_std) ** 2 * 2.0)
            result = np.where(np.isnan(result) | np.isinf(result), 1.0, result)
        return result.item() if np.ndim(result) == 0 else result