Skip to content

Xarray Accessors

The MONET Stats package provides first-class Xarray accessors to enable seamless integration with the Pangeo ecosystem. These accessors allow you to perform statistical analyses directly on xarray.DataArray and xarray.Dataset objects using the .monet_stats namespace.

DataArray Accessor

The MonetDataArrayAccessor (available via da.monet_stats) provides methods for common time-series and spatial analyses, as well as a comprehensive suite of verification metrics.

Verification Metrics

The accessor provides direct access to many common verification metrics. These methods typically take an obs DataArray as their first argument and an optional dim parameter.

Example usage:

rmse = mod_da.monet_stats.rmse(obs_da, dim="time")
pearson_r = mod_da.monet_stats.pearsonr(obs_da)

Available metrics include: - Error Metrics: mae, rmse, mb, ioa, crmse, mdnb, nmse, mnb, mne, nse - Correlation Metrics: pearsonr, r2, kge, ccc - Relative Metrics: nmb, fb

Verification Bundle

The verify method allows for efficient computation of multiple metrics at once, returning an xarray.Dataset.

metrics_ds = mod_da.monet_stats.verify(obs_da, dim="time")

The bundle includes: MAE, RMSE, MB, R (Pearson), IOA, NMB, MNB, MNE, NSE, and R2.

Full API Reference

Accessor for xarray.DataArray to provide MONET statistical methods.

Source code in src/monet_stats/accessor.py
 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
 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
113
114
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
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
229
230
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
287
288
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
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
390
391
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
459
460
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
527
528
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
605
606
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
683
684
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
753
754
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
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
@xr.register_dataarray_accessor("monet_stats")
class MonetDataArrayAccessor:
    """
    Accessor for xarray.DataArray to provide MONET statistical methods.
    """

    def __init__(self, xarray_obj: xr.DataArray) -> None:
        self._obj = xarray_obj

    def climatology(self, freq: str = "season", method: str = "mean", dim: str = "time") -> xr.DataArray:
        """
        Compute climatological statistics.

        Parameters
        ----------
        freq : str, optional
            Climatology frequency ('season', 'month', 'dayofyear', 'hour').
            Default is 'season'.
        method : str, optional
            Statistical method to apply ('mean', 'std', 'min', 'max', 'median').
            Default is 'mean'.
        dim : str, optional
            Dimension along which to compute climatology. Default is 'time'.

        Returns
        -------
        xarray.DataArray
            Climatological statistics.
        """
        return analysis.climatology(self._obj, freq=freq, method=method, dim=dim)

    def resample_data(self, freq: str = "MS", method: str = "mean", dim: str = "time", **kwargs: Any) -> xr.DataArray:
        """
        Resample data to a new temporal frequency.

        Parameters
        ----------
        freq : str, optional
            Resampling frequency (e.g., 'MS', 'W', 'D'). Default is 'MS'.
        method : str, optional
            Statistical method to apply ('mean', 'sum', 'min', 'max', 'std', 'median').
            Default is 'mean'.
        dim : str, optional
            Dimension along which to resample. Default is 'time'.
        **kwargs : Any
            Additional keyword arguments passed to the resample method.

        Returns
        -------
        xarray.DataArray
            Resampled data.
        """
        return analysis.resample_data(self._obj, freq=freq, method=method, dim=dim, **kwargs)

    def kz_filter(self, m: int, k: int, dim: str = "time") -> xr.DataArray:
        """
        Apply Kolmogorov-Zurbenko (KZ) filter.

        Parameters
        ----------
        m : int
            Window size for the moving average.
        k : int
            Number of iterations.
        dim : str, optional
            Dimension along which to apply the filter. Default is 'time'.

        Returns
        -------
        xarray.DataArray
            Filtered data.
        """
        return analysis.kz_filter(self._obj, m=m, k=k, dim=dim)

    def diurnal_cycle(self, method: str = "mean", dim: str = "time") -> xr.DataArray:
        """
        Compute the diurnal cycle (average hourly profile).

        Parameters
        ----------
        method : str, optional
            Statistical method to apply ('mean', 'median', 'std'). Default is 'mean'.
        dim : str, optional
            Dimension along which to compute the cycle. Default is 'time'.

        Returns
        -------
        xarray.DataArray
            Diurnal cycle.
        """
        return analysis.diurnal_cycle(self._obj, method=method, dim=dim)

    def rolling_mean_8h(self, dim: str = "time", min_periods: int = 6, center: bool = True) -> xr.DataArray:
        """
        Compute rolling 8-hour mean.

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute the mean. Default is 'time'.
        min_periods : int, optional
            Minimum number of observations in window. Default is 6.
        center : bool, optional
            If True, set the labels at the center of the window. Default is True.

        Returns
        -------
        xarray.DataArray
            Rolling 8-hour mean.
        """
        return analysis.rolling_mean_8h(self._obj, dim=dim, min_periods=min_periods, center=center)

    def rolling_mean_24h(self, dim: str = "time", min_periods: int = 18, center: bool = True) -> xr.DataArray:
        """
        Compute rolling 24-hour mean.

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute the mean. Default is 'time'.
        min_periods : int, optional
            Minimum number of observations in window. Default is 18.
        center : bool, optional
            If True, set the labels at the center of the window. Default is True.

        Returns
        -------
        xarray.DataArray
            Rolling 24-hour mean.
        """
        return analysis.rolling_mean_24h(self._obj, dim=dim, min_periods=min_periods, center=center)

    def mda1(self, dim: str = "time") -> xr.DataArray:
        """
        Compute Maximum Daily 1-hour Average (MDA1).

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute. Default is 'time'.

        Returns
        -------
        xarray.DataArray
            MDA1 values.
        """
        return analysis.mda1(self._obj, dim=dim)

    def mda8(self, dim: str = "time", min_periods: int = 6, center: bool = False) -> xr.DataArray:
        """
        Compute Maximum Daily 8-hour Average (MDA8).

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute. Default is 'time'.
        min_periods : int, optional
            Minimum number of observations for the 8-hour rolling mean. Default is 6.
        center : bool, optional
            Whether to center the 8-hour rolling window. Default is False.

        Returns
        -------
        xarray.DataArray
            MDA8 values.
        """
        return analysis.mda8(self._obj, dim=dim, min_periods=min_periods, center=center)

    def exceedance_count(self, threshold: float, dim: str = "time") -> xr.DataArray:
        """
        Count exceedances of a threshold.

        Parameters
        ----------
        threshold : float
            Value above which an exceedance is counted.
        dim : str, optional
            Dimension along which to count exceedances. Default is 'time'.

        Returns
        -------
        xarray.DataArray
            Number of exceedances.
        """
        return analysis.exceedance_count(self._obj, threshold=threshold, dim=dim)

    def percentile(self, q: Union[float, List[float], np.ndarray], dim: str = "time", **kwargs: Any) -> xr.DataArray:
        """
        Compute percentiles.

        Parameters
        ----------
        q : float or list of float
            Percentile(s) to compute (0-100).
        dim : str, optional
            Dimension over which to compute percentiles. Default is 'time'.
        **kwargs : Any
            Additional keyword arguments passed to xarray.quantile.

        Returns
        -------
        xarray.DataArray
            Computed percentiles.
        """
        return analysis.percentile(self._obj, q=q, dim=dim, **kwargs)

    def peak_timing(self, dim: str = "time") -> xr.DataArray:
        """
        Identify the coordinate value of the maximum.

        Parameters
        ----------
        dim : str, optional
            Dimension along which to find the peak. Default is 'time'.

        Returns
        -------
        xarray.DataArray
            Coordinate values where the maximum occurs.
        """
        return analysis.peak_timing(self._obj, dim=dim)

    def weighted_spatial_mean(
        self,
        lat_dim: str = "lat",
        lon_dim: str = "lon",
        weights: Optional[Union[xr.DataArray, np.ndarray]] = None,
    ) -> xr.DataArray:
        """
        Compute area-weighted spatial mean.

        Parameters
        ----------
        lat_dim : str, optional
            Name of the latitude dimension. Default is 'lat'.
        lon_dim : str, optional
            Name of the longitude dimension. Default is 'lon'.
        weights : xarray.DataArray or numpy.ndarray, optional
            Custom weights for the mean.

        Returns
        -------
        xarray.DataArray
            Area-weighted spatial mean.
        """
        return analysis.weighted_spatial_mean(self._obj, lat_dim=lat_dim, lon_dim=lon_dim, weights=weights)

    def fft_analysis(self, dim: str = "time", output: str = "psd") -> xr.DataArray:
        """
        Perform Fast Fourier Transform (FFT) analysis.

        Parameters
        ----------
        dim : str, optional
            Dimension along which to perform FFT. Default is 'time'.
        output : str, optional
            Type of output to return ('psd', 'magnitude', 'complex'). Default is 'psd'.

        Returns
        -------
        xarray.DataArray
            FFT results.
        """
        return analysis.fft_analysis(self._obj, dim=dim, output=output)

    def power_spectrum(
        self,
        dim: str = "time",
        fs: float = 1.0,
        window: str = "hann",
        nperseg: Optional[int] = None,
        **kwargs: Any,
    ) -> xr.DataArray:
        """
        Compute power spectrum using Welch's method.

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute the spectrum. Default is 'time'.
        fs : float, optional
            Sampling frequency. Default is 1.0.
        window : str, optional
            Desired window to use. Default is 'hann'.
        nperseg : int, optional
            Length of each segment.
        **kwargs : Any
            Additional keyword arguments passed to scipy.signal.welch.

        Returns
        -------
        xarray.DataArray
            Power spectral density.
        """
        return analysis.power_spectrum(self._obj, dim=dim, fs=fs, window=window, nperseg=nperseg, **kwargs)

    def seasonal_mean(self, dim: str = "time", weighted: bool = True) -> xr.DataArray:
        """
        Compute seasonal mean (DJF, MAM, JJA, SON).

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute the mean. Default is 'time'.
        weighted : bool, optional
            If True, weight by days in month. Default is True.

        Returns
        -------
        xarray.DataArray
            Seasonal means.
        """
        return analysis.seasonal_mean(self._obj, dim=dim, weighted=weighted)

    def monthly_climatology(self, dim: str = "time", method: str = "mean") -> xr.DataArray:
        """
        Compute monthly climatology.

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute the climatology. Default is 'time'.
        method : str, optional
            Statistical method to apply. Default is 'mean'.

        Returns
        -------
        xarray.DataArray
            Monthly climatology.
        """
        return analysis.monthly_climatology(self._obj, dim=dim, method=method)

    def anomalies(self, freq: str = "month", dim: str = "time") -> xr.DataArray:
        """
        Compute anomalies by subtracting the climatology.

        Parameters
        ----------
        freq : str, optional
            Climatology frequency ('season', 'month', 'dayofyear', 'hour').
            Default is 'month'.
        dim : str, optional
            Dimension along which to compute the anomalies. Default is 'time'.

        Returns
        -------
        xarray.DataArray
            Anomalies.
        """
        return analysis.anomalies(self._obj, freq=freq, dim=dim)

    def detrend(self, method: str = "linear", dim: str = "time") -> xr.DataArray:
        """
        Remove trend from data.

        Parameters
        ----------
        method : str, optional
            Detrending method ('linear', 'constant'). Default is 'linear'.
        dim : str, optional
            Dimension along which to detrend. Default is 'time'.

        Returns
        -------
        xarray.DataArray
            Detrended data.
        """
        return analysis.detrend(self._obj, method=method, dim=dim)

    def optimize(self, target_mb: float = 100.0) -> xr.DataArray:
        """
        Optimize performance by ensuring laziness and recommended chunks (Aero Protocol).

        Parameters
        ----------
        target_mb : float, optional
            Target size for each chunk in Megabytes. Default is 100.0.

        Returns
        -------
        xarray.DataArray
            Optimized DataArray.
        """
        from .utils_stats import _update_history

        if not performance._has_dask():
            return _update_history(self._obj, "Optimization skipped (Dask not installed)")

        # Ensure data is lazy
        res = performance.apply_lazy_threshold(self._obj, threshold_mb=0.1)
        # Always calculate and apply recommended chunks for the target size
        recommendation = performance.get_chunk_recommendation(res, target_mb=target_mb)
        res = res.chunk(recommendation)

        return _update_history(res, f"Optimized for performance (target={target_mb}MB)")

    def rechunk(self, chunks: Optional[dict] = None) -> xr.DataArray:
        """
        Apply new chunks to the DataArray (Aero Protocol provenance tracking).

        Parameters
        ----------
        chunks : dict, optional
            New chunk sizes. If None, uses optimal recommendations (~100MB).

        Returns
        -------
        xarray.DataArray
            Rechunked DataArray.
        """
        from .utils_stats import _update_history

        if not performance._has_dask():
            return _update_history(self._obj, "Rechunking skipped (Dask not installed)")

        if chunks is None:
            chunks = performance.get_chunk_recommendation(self._obj)

        res = self._obj.chunk(chunks)

        return _update_history(res, f"Rechunked with {chunks}")

    def taylor_statistics(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.Dataset:
        """
        Calculate components required for a Taylor diagram (Aero Protocol).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values (reference).
        dim : str or list of str, optional
            Dimension(s) along which to compute the statistics.

        Returns
        -------
        xarray.Dataset
            Dataset containing:
            - std_obs: Standard deviation of observations.
            - std_mod: Standard deviation of model predictions.
            - correlation: Pearson correlation coefficient.
        """
        obs, mod = xr.align(obs, self._obj, join="inner")
        from .utils_stats import _resolve_axis_to_dim, _update_history

        d = _resolve_axis_to_dim(obs, dim)

        res = xr.Dataset(
            {
                "std_obs": obs.std(dim=d),
                "std_mod": mod.std(dim=d),
                "correlation": correlation_metrics.pearsonr(obs, mod, axis=dim),
            }
        )
        return _update_history(res, "Taylor statistics components")

    def mae(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Mean Absolute Error (MAE).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Mean absolute error.
        """
        return error_metrics.MAE(obs, self._obj, axis=dim)

    def rmse(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Root Mean Square Error (RMSE).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Root mean square error.
        """
        return error_metrics.RMSE(obs, self._obj, axis=dim)

    def mb(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Mean Bias (MB).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Mean bias.
        """
        return error_metrics.MB(obs, self._obj, axis=dim)

    def ioa(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Index of Agreement (IOA).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Index of agreement.
        """
        return error_metrics.IOA(obs, self._obj, axis=dim)

    def crmse(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Centered Root Mean Square Error (CRMSE).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Centered root mean square error.
        """
        return error_metrics.CRMSE(obs, self._obj, axis=dim)

    def mdnb(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Median Bias (MdnB).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Median bias.
        """
        return error_metrics.MdnB(obs, self._obj, axis=dim)

    def nmse(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Normalized Mean Square Error (NMSE).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Normalized mean square error.
        """
        return error_metrics.NMSE(obs, self._obj, axis=dim)

    def pearsonr(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Pearson correlation coefficient.

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Pearson correlation coefficient.
        """
        return correlation_metrics.pearsonr(obs, self._obj, axis=dim)

    def r2(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Coefficient of Determination (R^2).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Coefficient of determination.
        """
        return correlation_metrics.R2(obs, self._obj, axis=dim)

    def kge(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Kling-Gupta Efficiency (KGE).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Kling-Gupta efficiency.
        """
        return correlation_metrics.KGE(obs, self._obj, axis=dim)

    def ccc(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Concordance Correlation Coefficient (CCC).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Concordance correlation coefficient.
        """
        return correlation_metrics.CCC(obs, self._obj, axis=dim)

    def nmb(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Normalized Mean Bias (NMB).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Normalized mean bias.
        """
        return relative_metrics.NMB(obs, self._obj, axis=dim)

    def fb(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Fractional Bias (FB).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Fractional bias.
        """
        return relative_metrics.FB(obs, self._obj, axis=dim)

    def mnb(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Mean Normalized Bias (MNB).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Mean normalized bias.
        """
        return error_metrics.MNB(obs, self._obj, axis=dim)

    def mne(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Mean Normalized Gross Error (MNE).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Mean normalized gross error.
        """
        return error_metrics.MNE(obs, self._obj, axis=dim)

    def nse(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
        """
        Compute Nash-Sutcliffe Efficiency (NSE).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metric.

        Returns
        -------
        xarray.DataArray
            Nash-Sutcliffe efficiency.
        """
        return efficiency_metrics.NSE(obs, self._obj, axis=dim)

    def verify(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.Dataset:
        """
        Calculate a bundle of common evaluation metrics (Aero Protocol).

        Parameters
        ----------
        obs : xarray.DataArray
            Observed values.
        dim : str or list of str, optional
            Dimension(s) along which to compute the metrics.

        Returns
        -------
        xarray.Dataset
            Dataset containing: MAE, RMSE, MB, R, IOA, NMB, MNB, MNE, NSE, and R2.
        """
        metrics = {
            "MAE": error_metrics.MAE(obs, self._obj, axis=dim),
            "RMSE": error_metrics.RMSE(obs, self._obj, axis=dim),
            "MB": error_metrics.MB(obs, self._obj, axis=dim),
            "R": correlation_metrics.pearsonr(obs, self._obj, axis=dim),
            "IOA": error_metrics.IOA(obs, self._obj, axis=dim),
            "NMB": relative_metrics.NMB(obs, self._obj, axis=dim),
            "MNB": error_metrics.MNB(obs, self._obj, axis=dim),
            "MNE": error_metrics.MNE(obs, self._obj, axis=dim),
            "NSE": efficiency_metrics.NSE(obs, self._obj, axis=dim),
            "R2": correlation_metrics.R2(obs, self._obj, axis=dim),
        }
        res = xr.Dataset(metrics)
        from .utils_stats import _update_history

        return _update_history(res, "Verification metrics bundle (verify)")

    def plot_spatial(
        self,
        method: str = "matplotlib",
        lat_dim: str = "lat",
        lon_dim: str = "lon",
        title: Optional[str] = None,
        cmap: str = "viridis",
        **kwargs: Any,
    ) -> Any:
        """
        Plot spatial data following the Aero Protocol's Two-Track Rule.

        Parameters
        ----------
        method : str, optional
            Plotting track: 'matplotlib' (Track A) or 'hvplot' (Track B).
            Default is 'matplotlib'.
        lat_dim : str, optional
            Latitude dimension name. Default is 'lat'.
        lon_dim : str, optional
            Longitude dimension name. Default is 'lon'.
        title : str, optional
            Plot title.
        cmap : str, optional
            Colormap. Default is 'viridis'.
        **kwargs : Any
            Additional keyword arguments.

        Returns
        -------
        Any
            The plot object.
        """
        from .visualize import plot_spatial

        return plot_spatial(
            self._obj,
            method=method,
            lat_dim=lat_dim,
            lon_dim=lon_dim,
            title=title,
            cmap=cmap,
            **kwargs,
        )

anomalies(freq='month', dim='time')

Compute anomalies by subtracting the climatology.

Parameters

freq : str, optional Climatology frequency ('season', 'month', 'dayofyear', 'hour'). Default is 'month'. dim : str, optional Dimension along which to compute the anomalies. Default is 'time'.

Returns

xarray.DataArray Anomalies.

Source code in src/monet_stats/accessor.py
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
def anomalies(self, freq: str = "month", dim: str = "time") -> xr.DataArray:
    """
    Compute anomalies by subtracting the climatology.

    Parameters
    ----------
    freq : str, optional
        Climatology frequency ('season', 'month', 'dayofyear', 'hour').
        Default is 'month'.
    dim : str, optional
        Dimension along which to compute the anomalies. Default is 'time'.

    Returns
    -------
    xarray.DataArray
        Anomalies.
    """
    return analysis.anomalies(self._obj, freq=freq, dim=dim)

ccc(obs, dim=None)

Compute Concordance Correlation Coefficient (CCC).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Concordance correlation coefficient.

Source code in src/monet_stats/accessor.py
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
def ccc(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Concordance Correlation Coefficient (CCC).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Concordance correlation coefficient.
    """
    return correlation_metrics.CCC(obs, self._obj, axis=dim)

climatology(freq='season', method='mean', dim='time')

Compute climatological statistics.

Parameters

freq : str, optional Climatology frequency ('season', 'month', 'dayofyear', 'hour'). Default is 'season'. method : str, optional Statistical method to apply ('mean', 'std', 'min', 'max', 'median'). Default is 'mean'. dim : str, optional Dimension along which to compute climatology. Default is 'time'.

Returns

xarray.DataArray Climatological statistics.

Source code in src/monet_stats/accessor.py
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
def climatology(self, freq: str = "season", method: str = "mean", dim: str = "time") -> xr.DataArray:
    """
    Compute climatological statistics.

    Parameters
    ----------
    freq : str, optional
        Climatology frequency ('season', 'month', 'dayofyear', 'hour').
        Default is 'season'.
    method : str, optional
        Statistical method to apply ('mean', 'std', 'min', 'max', 'median').
        Default is 'mean'.
    dim : str, optional
        Dimension along which to compute climatology. Default is 'time'.

    Returns
    -------
    xarray.DataArray
        Climatological statistics.
    """
    return analysis.climatology(self._obj, freq=freq, method=method, dim=dim)

crmse(obs, dim=None)

Compute Centered Root Mean Square Error (CRMSE).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Centered root mean square error.

Source code in src/monet_stats/accessor.py
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
def crmse(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Centered Root Mean Square Error (CRMSE).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Centered root mean square error.
    """
    return error_metrics.CRMSE(obs, self._obj, axis=dim)

detrend(method='linear', dim='time')

Remove trend from data.

Parameters

method : str, optional Detrending method ('linear', 'constant'). Default is 'linear'. dim : str, optional Dimension along which to detrend. Default is 'time'.

Returns

xarray.DataArray Detrended data.

Source code in src/monet_stats/accessor.py
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
def detrend(self, method: str = "linear", dim: str = "time") -> xr.DataArray:
    """
    Remove trend from data.

    Parameters
    ----------
    method : str, optional
        Detrending method ('linear', 'constant'). Default is 'linear'.
    dim : str, optional
        Dimension along which to detrend. Default is 'time'.

    Returns
    -------
    xarray.DataArray
        Detrended data.
    """
    return analysis.detrend(self._obj, method=method, dim=dim)

diurnal_cycle(method='mean', dim='time')

Compute the diurnal cycle (average hourly profile).

Parameters

method : str, optional Statistical method to apply ('mean', 'median', 'std'). Default is 'mean'. dim : str, optional Dimension along which to compute the cycle. Default is 'time'.

Returns

xarray.DataArray Diurnal cycle.

Source code in src/monet_stats/accessor.py
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
def diurnal_cycle(self, method: str = "mean", dim: str = "time") -> xr.DataArray:
    """
    Compute the diurnal cycle (average hourly profile).

    Parameters
    ----------
    method : str, optional
        Statistical method to apply ('mean', 'median', 'std'). Default is 'mean'.
    dim : str, optional
        Dimension along which to compute the cycle. Default is 'time'.

    Returns
    -------
    xarray.DataArray
        Diurnal cycle.
    """
    return analysis.diurnal_cycle(self._obj, method=method, dim=dim)

exceedance_count(threshold, dim='time')

Count exceedances of a threshold.

Parameters

threshold : float Value above which an exceedance is counted. dim : str, optional Dimension along which to count exceedances. Default is 'time'.

Returns

xarray.DataArray Number of exceedances.

Source code in src/monet_stats/accessor.py
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
def exceedance_count(self, threshold: float, dim: str = "time") -> xr.DataArray:
    """
    Count exceedances of a threshold.

    Parameters
    ----------
    threshold : float
        Value above which an exceedance is counted.
    dim : str, optional
        Dimension along which to count exceedances. Default is 'time'.

    Returns
    -------
    xarray.DataArray
        Number of exceedances.
    """
    return analysis.exceedance_count(self._obj, threshold=threshold, dim=dim)

fb(obs, dim=None)

Compute Fractional Bias (FB).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Fractional bias.

Source code in src/monet_stats/accessor.py
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
def fb(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Fractional Bias (FB).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Fractional bias.
    """
    return relative_metrics.FB(obs, self._obj, axis=dim)

fft_analysis(dim='time', output='psd')

Perform Fast Fourier Transform (FFT) analysis.

Parameters

dim : str, optional Dimension along which to perform FFT. Default is 'time'. output : str, optional Type of output to return ('psd', 'magnitude', 'complex'). Default is 'psd'.

Returns

xarray.DataArray FFT results.

Source code in src/monet_stats/accessor.py
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
def fft_analysis(self, dim: str = "time", output: str = "psd") -> xr.DataArray:
    """
    Perform Fast Fourier Transform (FFT) analysis.

    Parameters
    ----------
    dim : str, optional
        Dimension along which to perform FFT. Default is 'time'.
    output : str, optional
        Type of output to return ('psd', 'magnitude', 'complex'). Default is 'psd'.

    Returns
    -------
    xarray.DataArray
        FFT results.
    """
    return analysis.fft_analysis(self._obj, dim=dim, output=output)

ioa(obs, dim=None)

Compute Index of Agreement (IOA).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Index of agreement.

Source code in src/monet_stats/accessor.py
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
def ioa(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Index of Agreement (IOA).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Index of agreement.
    """
    return error_metrics.IOA(obs, self._obj, axis=dim)

kge(obs, dim=None)

Compute Kling-Gupta Efficiency (KGE).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Kling-Gupta efficiency.

Source code in src/monet_stats/accessor.py
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
def kge(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Kling-Gupta Efficiency (KGE).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Kling-Gupta efficiency.
    """
    return correlation_metrics.KGE(obs, self._obj, axis=dim)

kz_filter(m, k, dim='time')

Apply Kolmogorov-Zurbenko (KZ) filter.

Parameters

m : int Window size for the moving average. k : int Number of iterations. dim : str, optional Dimension along which to apply the filter. Default is 'time'.

Returns

xarray.DataArray Filtered data.

Source code in src/monet_stats/accessor.py
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
def kz_filter(self, m: int, k: int, dim: str = "time") -> xr.DataArray:
    """
    Apply Kolmogorov-Zurbenko (KZ) filter.

    Parameters
    ----------
    m : int
        Window size for the moving average.
    k : int
        Number of iterations.
    dim : str, optional
        Dimension along which to apply the filter. Default is 'time'.

    Returns
    -------
    xarray.DataArray
        Filtered data.
    """
    return analysis.kz_filter(self._obj, m=m, k=k, dim=dim)

mae(obs, dim=None)

Compute Mean Absolute Error (MAE).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Mean absolute error.

Source code in src/monet_stats/accessor.py
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
def mae(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Mean Absolute Error (MAE).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Mean absolute error.
    """
    return error_metrics.MAE(obs, self._obj, axis=dim)

mb(obs, dim=None)

Compute Mean Bias (MB).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Mean bias.

Source code in src/monet_stats/accessor.py
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
def mb(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Mean Bias (MB).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Mean bias.
    """
    return error_metrics.MB(obs, self._obj, axis=dim)

mda1(dim='time')

Compute Maximum Daily 1-hour Average (MDA1).

Parameters

dim : str, optional Dimension along which to compute. Default is 'time'.

Returns

xarray.DataArray MDA1 values.

Source code in src/monet_stats/accessor.py
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
def mda1(self, dim: str = "time") -> xr.DataArray:
    """
    Compute Maximum Daily 1-hour Average (MDA1).

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute. Default is 'time'.

    Returns
    -------
    xarray.DataArray
        MDA1 values.
    """
    return analysis.mda1(self._obj, dim=dim)

mda8(dim='time', min_periods=6, center=False)

Compute Maximum Daily 8-hour Average (MDA8).

Parameters

dim : str, optional Dimension along which to compute. Default is 'time'. min_periods : int, optional Minimum number of observations for the 8-hour rolling mean. Default is 6. center : bool, optional Whether to center the 8-hour rolling window. Default is False.

Returns

xarray.DataArray MDA8 values.

Source code in src/monet_stats/accessor.py
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
def mda8(self, dim: str = "time", min_periods: int = 6, center: bool = False) -> xr.DataArray:
    """
    Compute Maximum Daily 8-hour Average (MDA8).

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute. Default is 'time'.
    min_periods : int, optional
        Minimum number of observations for the 8-hour rolling mean. Default is 6.
    center : bool, optional
        Whether to center the 8-hour rolling window. Default is False.

    Returns
    -------
    xarray.DataArray
        MDA8 values.
    """
    return analysis.mda8(self._obj, dim=dim, min_periods=min_periods, center=center)

mdnb(obs, dim=None)

Compute Median Bias (MdnB).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Median bias.

Source code in src/monet_stats/accessor.py
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
def mdnb(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Median Bias (MdnB).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Median bias.
    """
    return error_metrics.MdnB(obs, self._obj, axis=dim)

mnb(obs, dim=None)

Compute Mean Normalized Bias (MNB).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Mean normalized bias.

Source code in src/monet_stats/accessor.py
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
def mnb(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Mean Normalized Bias (MNB).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Mean normalized bias.
    """
    return error_metrics.MNB(obs, self._obj, axis=dim)

mne(obs, dim=None)

Compute Mean Normalized Gross Error (MNE).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Mean normalized gross error.

Source code in src/monet_stats/accessor.py
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
def mne(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Mean Normalized Gross Error (MNE).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Mean normalized gross error.
    """
    return error_metrics.MNE(obs, self._obj, axis=dim)

monthly_climatology(dim='time', method='mean')

Compute monthly climatology.

Parameters

dim : str, optional Dimension along which to compute the climatology. Default is 'time'. method : str, optional Statistical method to apply. Default is 'mean'.

Returns

xarray.DataArray Monthly climatology.

Source code in src/monet_stats/accessor.py
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
def monthly_climatology(self, dim: str = "time", method: str = "mean") -> xr.DataArray:
    """
    Compute monthly climatology.

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute the climatology. Default is 'time'.
    method : str, optional
        Statistical method to apply. Default is 'mean'.

    Returns
    -------
    xarray.DataArray
        Monthly climatology.
    """
    return analysis.monthly_climatology(self._obj, dim=dim, method=method)

nmb(obs, dim=None)

Compute Normalized Mean Bias (NMB).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Normalized mean bias.

Source code in src/monet_stats/accessor.py
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
def nmb(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Normalized Mean Bias (NMB).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Normalized mean bias.
    """
    return relative_metrics.NMB(obs, self._obj, axis=dim)

nmse(obs, dim=None)

Compute Normalized Mean Square Error (NMSE).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Normalized mean square error.

Source code in src/monet_stats/accessor.py
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
def nmse(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Normalized Mean Square Error (NMSE).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Normalized mean square error.
    """
    return error_metrics.NMSE(obs, self._obj, axis=dim)

nse(obs, dim=None)

Compute Nash-Sutcliffe Efficiency (NSE).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Nash-Sutcliffe efficiency.

Source code in src/monet_stats/accessor.py
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
def nse(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Nash-Sutcliffe Efficiency (NSE).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Nash-Sutcliffe efficiency.
    """
    return efficiency_metrics.NSE(obs, self._obj, axis=dim)

optimize(target_mb=100.0)

Optimize performance by ensuring laziness and recommended chunks (Aero Protocol).

Parameters

target_mb : float, optional Target size for each chunk in Megabytes. Default is 100.0.

Returns

xarray.DataArray Optimized DataArray.

Source code in src/monet_stats/accessor.py
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
def optimize(self, target_mb: float = 100.0) -> xr.DataArray:
    """
    Optimize performance by ensuring laziness and recommended chunks (Aero Protocol).

    Parameters
    ----------
    target_mb : float, optional
        Target size for each chunk in Megabytes. Default is 100.0.

    Returns
    -------
    xarray.DataArray
        Optimized DataArray.
    """
    from .utils_stats import _update_history

    if not performance._has_dask():
        return _update_history(self._obj, "Optimization skipped (Dask not installed)")

    # Ensure data is lazy
    res = performance.apply_lazy_threshold(self._obj, threshold_mb=0.1)
    # Always calculate and apply recommended chunks for the target size
    recommendation = performance.get_chunk_recommendation(res, target_mb=target_mb)
    res = res.chunk(recommendation)

    return _update_history(res, f"Optimized for performance (target={target_mb}MB)")

peak_timing(dim='time')

Identify the coordinate value of the maximum.

Parameters

dim : str, optional Dimension along which to find the peak. Default is 'time'.

Returns

xarray.DataArray Coordinate values where the maximum occurs.

Source code in src/monet_stats/accessor.py
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
def peak_timing(self, dim: str = "time") -> xr.DataArray:
    """
    Identify the coordinate value of the maximum.

    Parameters
    ----------
    dim : str, optional
        Dimension along which to find the peak. Default is 'time'.

    Returns
    -------
    xarray.DataArray
        Coordinate values where the maximum occurs.
    """
    return analysis.peak_timing(self._obj, dim=dim)

pearsonr(obs, dim=None)

Compute Pearson correlation coefficient.

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Pearson correlation coefficient.

Source code in src/monet_stats/accessor.py
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
def pearsonr(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Pearson correlation coefficient.

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Pearson correlation coefficient.
    """
    return correlation_metrics.pearsonr(obs, self._obj, axis=dim)

percentile(q, dim='time', **kwargs)

Compute percentiles.

Parameters

q : float or list of float Percentile(s) to compute (0-100). dim : str, optional Dimension over which to compute percentiles. Default is 'time'. **kwargs : Any Additional keyword arguments passed to xarray.quantile.

Returns

xarray.DataArray Computed percentiles.

Source code in src/monet_stats/accessor.py
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
def percentile(self, q: Union[float, List[float], np.ndarray], dim: str = "time", **kwargs: Any) -> xr.DataArray:
    """
    Compute percentiles.

    Parameters
    ----------
    q : float or list of float
        Percentile(s) to compute (0-100).
    dim : str, optional
        Dimension over which to compute percentiles. Default is 'time'.
    **kwargs : Any
        Additional keyword arguments passed to xarray.quantile.

    Returns
    -------
    xarray.DataArray
        Computed percentiles.
    """
    return analysis.percentile(self._obj, q=q, dim=dim, **kwargs)

plot_spatial(method='matplotlib', lat_dim='lat', lon_dim='lon', title=None, cmap='viridis', **kwargs)

Plot spatial data following the Aero Protocol's Two-Track Rule.

Parameters

method : str, optional Plotting track: 'matplotlib' (Track A) or 'hvplot' (Track B). Default is 'matplotlib'. lat_dim : str, optional Latitude dimension name. Default is 'lat'. lon_dim : str, optional Longitude dimension name. Default is 'lon'. title : str, optional Plot title. cmap : str, optional Colormap. Default is 'viridis'. **kwargs : Any Additional keyword arguments.

Returns

Any The plot object.

Source code in src/monet_stats/accessor.py
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
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
def plot_spatial(
    self,
    method: str = "matplotlib",
    lat_dim: str = "lat",
    lon_dim: str = "lon",
    title: Optional[str] = None,
    cmap: str = "viridis",
    **kwargs: Any,
) -> Any:
    """
    Plot spatial data following the Aero Protocol's Two-Track Rule.

    Parameters
    ----------
    method : str, optional
        Plotting track: 'matplotlib' (Track A) or 'hvplot' (Track B).
        Default is 'matplotlib'.
    lat_dim : str, optional
        Latitude dimension name. Default is 'lat'.
    lon_dim : str, optional
        Longitude dimension name. Default is 'lon'.
    title : str, optional
        Plot title.
    cmap : str, optional
        Colormap. Default is 'viridis'.
    **kwargs : Any
        Additional keyword arguments.

    Returns
    -------
    Any
        The plot object.
    """
    from .visualize import plot_spatial

    return plot_spatial(
        self._obj,
        method=method,
        lat_dim=lat_dim,
        lon_dim=lon_dim,
        title=title,
        cmap=cmap,
        **kwargs,
    )

power_spectrum(dim='time', fs=1.0, window='hann', nperseg=None, **kwargs)

Compute power spectrum using Welch's method.

Parameters

dim : str, optional Dimension along which to compute the spectrum. Default is 'time'. fs : float, optional Sampling frequency. Default is 1.0. window : str, optional Desired window to use. Default is 'hann'. nperseg : int, optional Length of each segment. **kwargs : Any Additional keyword arguments passed to scipy.signal.welch.

Returns

xarray.DataArray Power spectral density.

Source code in src/monet_stats/accessor.py
285
286
287
288
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
def power_spectrum(
    self,
    dim: str = "time",
    fs: float = 1.0,
    window: str = "hann",
    nperseg: Optional[int] = None,
    **kwargs: Any,
) -> xr.DataArray:
    """
    Compute power spectrum using Welch's method.

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute the spectrum. Default is 'time'.
    fs : float, optional
        Sampling frequency. Default is 1.0.
    window : str, optional
        Desired window to use. Default is 'hann'.
    nperseg : int, optional
        Length of each segment.
    **kwargs : Any
        Additional keyword arguments passed to scipy.signal.welch.

    Returns
    -------
    xarray.DataArray
        Power spectral density.
    """
    return analysis.power_spectrum(self._obj, dim=dim, fs=fs, window=window, nperseg=nperseg, **kwargs)

r2(obs, dim=None)

Compute Coefficient of Determination (R^2).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Coefficient of determination.

Source code in src/monet_stats/accessor.py
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
def r2(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Coefficient of Determination (R^2).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Coefficient of determination.
    """
    return correlation_metrics.R2(obs, self._obj, axis=dim)

rechunk(chunks=None)

Apply new chunks to the DataArray (Aero Protocol provenance tracking).

Parameters

chunks : dict, optional New chunk sizes. If None, uses optimal recommendations (~100MB).

Returns

xarray.DataArray Rechunked DataArray.

Source code in src/monet_stats/accessor.py
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
def rechunk(self, chunks: Optional[dict] = None) -> xr.DataArray:
    """
    Apply new chunks to the DataArray (Aero Protocol provenance tracking).

    Parameters
    ----------
    chunks : dict, optional
        New chunk sizes. If None, uses optimal recommendations (~100MB).

    Returns
    -------
    xarray.DataArray
        Rechunked DataArray.
    """
    from .utils_stats import _update_history

    if not performance._has_dask():
        return _update_history(self._obj, "Rechunking skipped (Dask not installed)")

    if chunks is None:
        chunks = performance.get_chunk_recommendation(self._obj)

    res = self._obj.chunk(chunks)

    return _update_history(res, f"Rechunked with {chunks}")

resample_data(freq='MS', method='mean', dim='time', **kwargs)

Resample data to a new temporal frequency.

Parameters

freq : str, optional Resampling frequency (e.g., 'MS', 'W', 'D'). Default is 'MS'. method : str, optional Statistical method to apply ('mean', 'sum', 'min', 'max', 'std', 'median'). Default is 'mean'. dim : str, optional Dimension along which to resample. Default is 'time'. **kwargs : Any Additional keyword arguments passed to the resample method.

Returns

xarray.DataArray Resampled data.

Source code in src/monet_stats/accessor.py
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
def resample_data(self, freq: str = "MS", method: str = "mean", dim: str = "time", **kwargs: Any) -> xr.DataArray:
    """
    Resample data to a new temporal frequency.

    Parameters
    ----------
    freq : str, optional
        Resampling frequency (e.g., 'MS', 'W', 'D'). Default is 'MS'.
    method : str, optional
        Statistical method to apply ('mean', 'sum', 'min', 'max', 'std', 'median').
        Default is 'mean'.
    dim : str, optional
        Dimension along which to resample. Default is 'time'.
    **kwargs : Any
        Additional keyword arguments passed to the resample method.

    Returns
    -------
    xarray.DataArray
        Resampled data.
    """
    return analysis.resample_data(self._obj, freq=freq, method=method, dim=dim, **kwargs)

rmse(obs, dim=None)

Compute Root Mean Square Error (RMSE).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metric.

Returns

xarray.DataArray Root mean square error.

Source code in src/monet_stats/accessor.py
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
def rmse(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.DataArray:
    """
    Compute Root Mean Square Error (RMSE).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metric.

    Returns
    -------
    xarray.DataArray
        Root mean square error.
    """
    return error_metrics.RMSE(obs, self._obj, axis=dim)

rolling_mean_24h(dim='time', min_periods=18, center=True)

Compute rolling 24-hour mean.

Parameters

dim : str, optional Dimension along which to compute the mean. Default is 'time'. min_periods : int, optional Minimum number of observations in window. Default is 18. center : bool, optional If True, set the labels at the center of the window. Default is True.

Returns

xarray.DataArray Rolling 24-hour mean.

Source code in src/monet_stats/accessor.py
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
def rolling_mean_24h(self, dim: str = "time", min_periods: int = 18, center: bool = True) -> xr.DataArray:
    """
    Compute rolling 24-hour mean.

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute the mean. Default is 'time'.
    min_periods : int, optional
        Minimum number of observations in window. Default is 18.
    center : bool, optional
        If True, set the labels at the center of the window. Default is True.

    Returns
    -------
    xarray.DataArray
        Rolling 24-hour mean.
    """
    return analysis.rolling_mean_24h(self._obj, dim=dim, min_periods=min_periods, center=center)

rolling_mean_8h(dim='time', min_periods=6, center=True)

Compute rolling 8-hour mean.

Parameters

dim : str, optional Dimension along which to compute the mean. Default is 'time'. min_periods : int, optional Minimum number of observations in window. Default is 6. center : bool, optional If True, set the labels at the center of the window. Default is True.

Returns

xarray.DataArray Rolling 8-hour mean.

Source code in src/monet_stats/accessor.py
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
def rolling_mean_8h(self, dim: str = "time", min_periods: int = 6, center: bool = True) -> xr.DataArray:
    """
    Compute rolling 8-hour mean.

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute the mean. Default is 'time'.
    min_periods : int, optional
        Minimum number of observations in window. Default is 6.
    center : bool, optional
        If True, set the labels at the center of the window. Default is True.

    Returns
    -------
    xarray.DataArray
        Rolling 8-hour mean.
    """
    return analysis.rolling_mean_8h(self._obj, dim=dim, min_periods=min_periods, center=center)

seasonal_mean(dim='time', weighted=True)

Compute seasonal mean (DJF, MAM, JJA, SON).

Parameters

dim : str, optional Dimension along which to compute the mean. Default is 'time'. weighted : bool, optional If True, weight by days in month. Default is True.

Returns

xarray.DataArray Seasonal means.

Source code in src/monet_stats/accessor.py
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
def seasonal_mean(self, dim: str = "time", weighted: bool = True) -> xr.DataArray:
    """
    Compute seasonal mean (DJF, MAM, JJA, SON).

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute the mean. Default is 'time'.
    weighted : bool, optional
        If True, weight by days in month. Default is True.

    Returns
    -------
    xarray.DataArray
        Seasonal means.
    """
    return analysis.seasonal_mean(self._obj, dim=dim, weighted=weighted)

taylor_statistics(obs, dim=None)

Calculate components required for a Taylor diagram (Aero Protocol).

Parameters

obs : xarray.DataArray Observed values (reference). dim : str or list of str, optional Dimension(s) along which to compute the statistics.

Returns

xarray.Dataset Dataset containing: - std_obs: Standard deviation of observations. - std_mod: Standard deviation of model predictions. - correlation: Pearson correlation coefficient.

Source code in src/monet_stats/accessor.py
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
def taylor_statistics(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.Dataset:
    """
    Calculate components required for a Taylor diagram (Aero Protocol).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values (reference).
    dim : str or list of str, optional
        Dimension(s) along which to compute the statistics.

    Returns
    -------
    xarray.Dataset
        Dataset containing:
        - std_obs: Standard deviation of observations.
        - std_mod: Standard deviation of model predictions.
        - correlation: Pearson correlation coefficient.
    """
    obs, mod = xr.align(obs, self._obj, join="inner")
    from .utils_stats import _resolve_axis_to_dim, _update_history

    d = _resolve_axis_to_dim(obs, dim)

    res = xr.Dataset(
        {
            "std_obs": obs.std(dim=d),
            "std_mod": mod.std(dim=d),
            "correlation": correlation_metrics.pearsonr(obs, mod, axis=dim),
        }
    )
    return _update_history(res, "Taylor statistics components")

verify(obs, dim=None)

Calculate a bundle of common evaluation metrics (Aero Protocol).

Parameters

obs : xarray.DataArray Observed values. dim : str or list of str, optional Dimension(s) along which to compute the metrics.

Returns

xarray.Dataset Dataset containing: MAE, RMSE, MB, R, IOA, NMB, MNB, MNE, NSE, and R2.

Source code in src/monet_stats/accessor.py
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
def verify(self, obs: xr.DataArray, dim: Optional[Union[str, List[str]]] = None) -> xr.Dataset:
    """
    Calculate a bundle of common evaluation metrics (Aero Protocol).

    Parameters
    ----------
    obs : xarray.DataArray
        Observed values.
    dim : str or list of str, optional
        Dimension(s) along which to compute the metrics.

    Returns
    -------
    xarray.Dataset
        Dataset containing: MAE, RMSE, MB, R, IOA, NMB, MNB, MNE, NSE, and R2.
    """
    metrics = {
        "MAE": error_metrics.MAE(obs, self._obj, axis=dim),
        "RMSE": error_metrics.RMSE(obs, self._obj, axis=dim),
        "MB": error_metrics.MB(obs, self._obj, axis=dim),
        "R": correlation_metrics.pearsonr(obs, self._obj, axis=dim),
        "IOA": error_metrics.IOA(obs, self._obj, axis=dim),
        "NMB": relative_metrics.NMB(obs, self._obj, axis=dim),
        "MNB": error_metrics.MNB(obs, self._obj, axis=dim),
        "MNE": error_metrics.MNE(obs, self._obj, axis=dim),
        "NSE": efficiency_metrics.NSE(obs, self._obj, axis=dim),
        "R2": correlation_metrics.R2(obs, self._obj, axis=dim),
    }
    res = xr.Dataset(metrics)
    from .utils_stats import _update_history

    return _update_history(res, "Verification metrics bundle (verify)")

weighted_spatial_mean(lat_dim='lat', lon_dim='lon', weights=None)

Compute area-weighted spatial mean.

Parameters

lat_dim : str, optional Name of the latitude dimension. Default is 'lat'. lon_dim : str, optional Name of the longitude dimension. Default is 'lon'. weights : xarray.DataArray or numpy.ndarray, optional Custom weights for the mean.

Returns

xarray.DataArray Area-weighted spatial mean.

Source code in src/monet_stats/accessor.py
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
def weighted_spatial_mean(
    self,
    lat_dim: str = "lat",
    lon_dim: str = "lon",
    weights: Optional[Union[xr.DataArray, np.ndarray]] = None,
) -> xr.DataArray:
    """
    Compute area-weighted spatial mean.

    Parameters
    ----------
    lat_dim : str, optional
        Name of the latitude dimension. Default is 'lat'.
    lon_dim : str, optional
        Name of the longitude dimension. Default is 'lon'.
    weights : xarray.DataArray or numpy.ndarray, optional
        Custom weights for the mean.

    Returns
    -------
    xarray.DataArray
        Area-weighted spatial mean.
    """
    return analysis.weighted_spatial_mean(self._obj, lat_dim=lat_dim, lon_dim=lon_dim, weights=weights)

Dataset Accessor

The MonetDatasetAccessor (available via ds.monet_stats) provides methods for calculating summary statistics and performing analyses across multiple variables in a dataset.

Accessor for xarray.Dataset to provide MONET statistical methods.

Source code in src/monet_stats/accessor.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
 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
 919
 920
 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
1008
1009
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
1079
1080
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
1176
1177
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
@xr.register_dataset_accessor("monet_stats")
class MonetDatasetAccessor:
    """
    Accessor for xarray.Dataset to provide MONET statistical methods.
    """

    def __init__(self, xarray_obj: xr.Dataset) -> None:
        self._obj = xarray_obj

    def stats(
        self,
        obs_name: str = "Obs",
        mod_name: str = "Mod",
        threshold: float = 0.0,
        minval: Optional[float] = None,
        maxval: Optional[float] = None,
    ) -> dict:
        """
        Calculate summary statistics for observations and model results.

        Parameters
        ----------
        obs_name : str, optional
            Name of observation variable. Default is 'Obs'.
        mod_name : str, optional
            Name of model variable. Default is 'Mod'.
        threshold : float, optional
            Threshold for contingency scores. Default is 0.0.
        minval : float, optional
            Minimum value for filtering.
        maxval : float, optional
            Maximum value for filtering.

        Returns
        -------
        dict
            Dictionary of calculated statistics.
        """
        from . import stats

        return stats(
            self._obj,
            obs_name=obs_name,
            mod_name=mod_name,
            threshold=threshold,
            minval=minval,
            maxval=maxval,
        )

    def climatology(self, freq: str = "season", method: str = "mean", dim: str = "time") -> xr.Dataset:
        """
        Compute climatological statistics.

        Parameters
        ----------
        freq : str, optional
            Climatology frequency. Default is 'season'.
        method : str, optional
            Statistical method. Default is 'mean'.
        dim : str, optional
            Dimension along which to compute. Default is 'time'.

        Returns
        -------
        xarray.Dataset
            Climatological statistics.
        """
        return analysis.climatology(self._obj, freq=freq, method=method, dim=dim)

    def resample_data(self, freq: str = "MS", method: str = "mean", dim: str = "time", **kwargs: Any) -> xr.Dataset:
        """
        Resample data to a new temporal frequency.

        Parameters
        ----------
        freq : str, optional
            Resampling frequency. Default is 'MS'.
        method : str, optional
            Statistical method. Default is 'mean'.
        dim : str, optional
            Dimension along which to resample. Default is 'time'.
        **kwargs : Any
            Additional keyword arguments.

        Returns
        -------
        xarray.Dataset
            Resampled data.
        """
        return analysis.resample_data(self._obj, freq=freq, method=method, dim=dim, **kwargs)

    def kz_filter(self, m: int, k: int, dim: str = "time") -> xr.Dataset:
        """
        Apply Kolmogorov-Zurbenko (KZ) filter.

        Parameters
        ----------
        m : int
            Window size.
        k : int
            Number of iterations.
        dim : str, optional
            Dimension along which to apply. Default is 'time'.

        Returns
        -------
        xarray.Dataset
            Filtered data.
        """
        return analysis.kz_filter(self._obj, m=m, k=k, dim=dim)

    def diurnal_cycle(self, method: str = "mean", dim: str = "time") -> xr.Dataset:
        """
        Compute the diurnal cycle.

        Parameters
        ----------
        method : str, optional
            Statistical method. Default is 'mean'.
        dim : str, optional
            Dimension along which to compute. Default is 'time'.

        Returns
        -------
        xarray.Dataset
            Diurnal cycle.
        """
        return analysis.diurnal_cycle(self._obj, method=method, dim=dim)

    def rolling_mean_8h(self, dim: str = "time", min_periods: int = 6, center: bool = True) -> xr.Dataset:
        """
        Compute rolling 8-hour mean.

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute. Default is 'time'.
        min_periods : int, optional
            Minimum number of observations. Default is 6.
        center : bool, optional
            If True, center the labels. Default is True.

        Returns
        -------
        xarray.Dataset
            Rolling 8-hour mean.
        """
        return analysis.rolling_mean_8h(self._obj, dim=dim, min_periods=min_periods, center=center)

    def rolling_mean_24h(self, dim: str = "time", min_periods: int = 18, center: bool = True) -> xr.Dataset:
        """
        Compute rolling 24-hour mean.

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute. Default is 'time'.
        min_periods : int, optional
            Minimum number of observations. Default is 18.
        center : bool, optional
            If True, center the labels. Default is True.

        Returns
        -------
        xarray.Dataset
            Rolling 24-hour mean.
        """
        return analysis.rolling_mean_24h(self._obj, dim=dim, min_periods=min_periods, center=center)

    def mda1(self, dim: str = "time") -> xr.Dataset:
        """
        Compute Maximum Daily 1-hour Average (MDA1).

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute. Default is 'time'.

        Returns
        -------
        xarray.Dataset
            MDA1 values.
        """
        return analysis.mda1(self._obj, dim=dim)

    def mda8(self, dim: str = "time", min_periods: int = 6, center: bool = False) -> xr.Dataset:
        """
        Compute Maximum Daily 8-hour Average (MDA8).

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute. Default is 'time'.
        min_periods : int, optional
            Minimum number of observations. Default is 6.
        center : bool, optional
            Whether to center the window. Default is False.

        Returns
        -------
        xarray.Dataset
            MDA8 values.
        """
        return analysis.mda8(self._obj, dim=dim, min_periods=min_periods, center=center)

    def exceedance_count(self, threshold: float, dim: str = "time") -> xr.Dataset:
        """
        Count exceedances of a threshold.

        Parameters
        ----------
        threshold : float
            Threshold value.
        dim : str, optional
            Dimension along which to count. Default is 'time'.

        Returns
        -------
        xarray.Dataset
            Number of exceedances.
        """
        return analysis.exceedance_count(self._obj, threshold=threshold, dim=dim)

    def percentile(self, q: Union[float, List[float], np.ndarray], dim: str = "time", **kwargs: Any) -> xr.Dataset:
        """
        Compute percentiles.

        Parameters
        ----------
        q : float or list of float
            Percentile(s) (0-100).
        dim : str, optional
            Dimension over which to compute. Default is 'time'.
        **kwargs : Any
            Additional keyword arguments.

        Returns
        -------
        xarray.Dataset
            Computed percentiles.
        """
        return analysis.percentile(self._obj, q=q, dim=dim, **kwargs)

    def peak_timing(self, dim: str = "time") -> xr.Dataset:
        """
        Identify the coordinate value of the maximum.

        Parameters
        ----------
        dim : str, optional
            Dimension along which to find peak. Default is 'time'.

        Returns
        -------
        xarray.Dataset
            Coordinate values.
        """
        return analysis.peak_timing(self._obj, dim=dim)

    def weighted_spatial_mean(
        self,
        lat_dim: str = "lat",
        lon_dim: str = "lon",
        weights: Optional[Union[xr.DataArray, np.ndarray]] = None,
    ) -> xr.Dataset:
        """
        Compute area-weighted spatial mean.

        Parameters
        ----------
        lat_dim : str, optional
            Latitude dimension name. Default is 'lat'.
        lon_dim : str, optional
            Longitude dimension name. Default is 'lon'.
        weights : xarray.DataArray or numpy.ndarray, optional
            Custom weights.

        Returns
        -------
        xarray.Dataset
            Area-weighted spatial mean.
        """
        return analysis.weighted_spatial_mean(self._obj, lat_dim=lat_dim, lon_dim=lon_dim, weights=weights)

    def seasonal_mean(self, dim: str = "time", weighted: bool = True) -> xr.Dataset:
        """
        Compute seasonal mean (DJF, MAM, JJA, SON).

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute. Default is 'time'.
        weighted : bool, optional
            Weight by days in month. Default is True.

        Returns
        -------
        xarray.Dataset
            Seasonal means.
        """
        return analysis.seasonal_mean(self._obj, dim=dim, weighted=weighted)

    def monthly_climatology(self, dim: str = "time", method: str = "mean") -> xr.Dataset:
        """
        Compute monthly climatology.

        Parameters
        ----------
        dim : str, optional
            Dimension along which to compute. Default is 'time'.
        method : str, optional
            Statistical method. Default is 'mean'.

        Returns
        -------
        xarray.Dataset
            Monthly climatology.
        """
        return analysis.monthly_climatology(self._obj, dim=dim, method=method)

    def anomalies(self, freq: str = "month", dim: str = "time") -> xr.Dataset:
        """
        Compute anomalies by subtracting the climatology.

        Parameters
        ----------
        freq : str, optional
            Climatology frequency ('season', 'month', 'dayofyear', 'hour').
            Default is 'month'.
        dim : str, optional
            Dimension along which to compute the anomalies. Default is 'time'.

        Returns
        -------
        xarray.Dataset
            Anomalies.
        """
        return analysis.anomalies(self._obj, freq=freq, dim=dim)

    def detrend(self, method: str = "linear", dim: str = "time") -> xr.Dataset:
        """
        Remove trend from data.

        Parameters
        ----------
        method : str, optional
            Detrending method ('linear', 'constant'). Default is 'linear'.
        dim : str, optional
            Dimension along which to detrend. Default is 'time'.

        Returns
        -------
        xarray.Dataset
            Detrended data.
        """
        return analysis.detrend(self._obj, method=method, dim=dim)

    def optimize(self, target_mb: float = 100.0) -> xr.Dataset:
        """
        Optimize performance by ensuring laziness and recommended chunks (Aero Protocol).

        Parameters
        ----------
        target_mb : float, optional
            Target size for each chunk in Megabytes. Default is 100.0.

        Returns
        -------
        xarray.Dataset
            Optimized Dataset.
        """
        from .utils_stats import _update_history

        if not performance._has_dask():
            return _update_history(self._obj, "Optimization skipped (Dask not installed)")

        # Ensure data is lazy
        res = performance.apply_lazy_threshold(self._obj, threshold_mb=0.1)
        # Always calculate and apply recommended chunks for the target size
        recommendation = performance.get_chunk_recommendation(res, target_mb=target_mb)
        res = res.chunk(recommendation)

        return _update_history(res, f"Optimized for performance (target={target_mb}MB)")

    def rechunk(self, chunks: Optional[dict] = None) -> xr.Dataset:
        """
        Apply new chunks to the Dataset (Aero Protocol provenance tracking).

        Parameters
        ----------
        chunks : dict, optional
            New chunk sizes. If None, uses optimal recommendations (~100MB).

        Returns
        -------
        xarray.Dataset
            Rechunked Dataset.
        """
        from .utils_stats import _update_history

        if not performance._has_dask():
            return _update_history(self._obj, "Rechunking skipped (Dask not installed)")

        if chunks is None:
            chunks = performance.get_chunk_recommendation(self._obj)

        res = self._obj.chunk(chunks)

        return _update_history(res, f"Rechunked with {chunks}")

anomalies(freq='month', dim='time')

Compute anomalies by subtracting the climatology.

Parameters

freq : str, optional Climatology frequency ('season', 'month', 'dayofyear', 'hour'). Default is 'month'. dim : str, optional Dimension along which to compute the anomalies. Default is 'time'.

Returns

xarray.Dataset Anomalies.

Source code in src/monet_stats/accessor.py
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
def anomalies(self, freq: str = "month", dim: str = "time") -> xr.Dataset:
    """
    Compute anomalies by subtracting the climatology.

    Parameters
    ----------
    freq : str, optional
        Climatology frequency ('season', 'month', 'dayofyear', 'hour').
        Default is 'month'.
    dim : str, optional
        Dimension along which to compute the anomalies. Default is 'time'.

    Returns
    -------
    xarray.Dataset
        Anomalies.
    """
    return analysis.anomalies(self._obj, freq=freq, dim=dim)

climatology(freq='season', method='mean', dim='time')

Compute climatological statistics.

Parameters

freq : str, optional Climatology frequency. Default is 'season'. method : str, optional Statistical method. Default is 'mean'. dim : str, optional Dimension along which to compute. Default is 'time'.

Returns

xarray.Dataset Climatological statistics.

Source code in src/monet_stats/accessor.py
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
def climatology(self, freq: str = "season", method: str = "mean", dim: str = "time") -> xr.Dataset:
    """
    Compute climatological statistics.

    Parameters
    ----------
    freq : str, optional
        Climatology frequency. Default is 'season'.
    method : str, optional
        Statistical method. Default is 'mean'.
    dim : str, optional
        Dimension along which to compute. Default is 'time'.

    Returns
    -------
    xarray.Dataset
        Climatological statistics.
    """
    return analysis.climatology(self._obj, freq=freq, method=method, dim=dim)

detrend(method='linear', dim='time')

Remove trend from data.

Parameters

method : str, optional Detrending method ('linear', 'constant'). Default is 'linear'. dim : str, optional Dimension along which to detrend. Default is 'time'.

Returns

xarray.Dataset Detrended data.

Source code in src/monet_stats/accessor.py
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
def detrend(self, method: str = "linear", dim: str = "time") -> xr.Dataset:
    """
    Remove trend from data.

    Parameters
    ----------
    method : str, optional
        Detrending method ('linear', 'constant'). Default is 'linear'.
    dim : str, optional
        Dimension along which to detrend. Default is 'time'.

    Returns
    -------
    xarray.Dataset
        Detrended data.
    """
    return analysis.detrend(self._obj, method=method, dim=dim)

diurnal_cycle(method='mean', dim='time')

Compute the diurnal cycle.

Parameters

method : str, optional Statistical method. Default is 'mean'. dim : str, optional Dimension along which to compute. Default is 'time'.

Returns

xarray.Dataset Diurnal cycle.

Source code in src/monet_stats/accessor.py
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
def diurnal_cycle(self, method: str = "mean", dim: str = "time") -> xr.Dataset:
    """
    Compute the diurnal cycle.

    Parameters
    ----------
    method : str, optional
        Statistical method. Default is 'mean'.
    dim : str, optional
        Dimension along which to compute. Default is 'time'.

    Returns
    -------
    xarray.Dataset
        Diurnal cycle.
    """
    return analysis.diurnal_cycle(self._obj, method=method, dim=dim)

exceedance_count(threshold, dim='time')

Count exceedances of a threshold.

Parameters

threshold : float Threshold value. dim : str, optional Dimension along which to count. Default is 'time'.

Returns

xarray.Dataset Number of exceedances.

Source code in src/monet_stats/accessor.py
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
def exceedance_count(self, threshold: float, dim: str = "time") -> xr.Dataset:
    """
    Count exceedances of a threshold.

    Parameters
    ----------
    threshold : float
        Threshold value.
    dim : str, optional
        Dimension along which to count. Default is 'time'.

    Returns
    -------
    xarray.Dataset
        Number of exceedances.
    """
    return analysis.exceedance_count(self._obj, threshold=threshold, dim=dim)

kz_filter(m, k, dim='time')

Apply Kolmogorov-Zurbenko (KZ) filter.

Parameters

m : int Window size. k : int Number of iterations. dim : str, optional Dimension along which to apply. Default is 'time'.

Returns

xarray.Dataset Filtered data.

Source code in src/monet_stats/accessor.py
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
def kz_filter(self, m: int, k: int, dim: str = "time") -> xr.Dataset:
    """
    Apply Kolmogorov-Zurbenko (KZ) filter.

    Parameters
    ----------
    m : int
        Window size.
    k : int
        Number of iterations.
    dim : str, optional
        Dimension along which to apply. Default is 'time'.

    Returns
    -------
    xarray.Dataset
        Filtered data.
    """
    return analysis.kz_filter(self._obj, m=m, k=k, dim=dim)

mda1(dim='time')

Compute Maximum Daily 1-hour Average (MDA1).

Parameters

dim : str, optional Dimension along which to compute. Default is 'time'.

Returns

xarray.Dataset MDA1 values.

Source code in src/monet_stats/accessor.py
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
def mda1(self, dim: str = "time") -> xr.Dataset:
    """
    Compute Maximum Daily 1-hour Average (MDA1).

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute. Default is 'time'.

    Returns
    -------
    xarray.Dataset
        MDA1 values.
    """
    return analysis.mda1(self._obj, dim=dim)

mda8(dim='time', min_periods=6, center=False)

Compute Maximum Daily 8-hour Average (MDA8).

Parameters

dim : str, optional Dimension along which to compute. Default is 'time'. min_periods : int, optional Minimum number of observations. Default is 6. center : bool, optional Whether to center the window. Default is False.

Returns

xarray.Dataset MDA8 values.

Source code in src/monet_stats/accessor.py
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
def mda8(self, dim: str = "time", min_periods: int = 6, center: bool = False) -> xr.Dataset:
    """
    Compute Maximum Daily 8-hour Average (MDA8).

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute. Default is 'time'.
    min_periods : int, optional
        Minimum number of observations. Default is 6.
    center : bool, optional
        Whether to center the window. Default is False.

    Returns
    -------
    xarray.Dataset
        MDA8 values.
    """
    return analysis.mda8(self._obj, dim=dim, min_periods=min_periods, center=center)

monthly_climatology(dim='time', method='mean')

Compute monthly climatology.

Parameters

dim : str, optional Dimension along which to compute. Default is 'time'. method : str, optional Statistical method. Default is 'mean'.

Returns

xarray.Dataset Monthly climatology.

Source code in src/monet_stats/accessor.py
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
def monthly_climatology(self, dim: str = "time", method: str = "mean") -> xr.Dataset:
    """
    Compute monthly climatology.

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute. Default is 'time'.
    method : str, optional
        Statistical method. Default is 'mean'.

    Returns
    -------
    xarray.Dataset
        Monthly climatology.
    """
    return analysis.monthly_climatology(self._obj, dim=dim, method=method)

optimize(target_mb=100.0)

Optimize performance by ensuring laziness and recommended chunks (Aero Protocol).

Parameters

target_mb : float, optional Target size for each chunk in Megabytes. Default is 100.0.

Returns

xarray.Dataset Optimized Dataset.

Source code in src/monet_stats/accessor.py
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
def optimize(self, target_mb: float = 100.0) -> xr.Dataset:
    """
    Optimize performance by ensuring laziness and recommended chunks (Aero Protocol).

    Parameters
    ----------
    target_mb : float, optional
        Target size for each chunk in Megabytes. Default is 100.0.

    Returns
    -------
    xarray.Dataset
        Optimized Dataset.
    """
    from .utils_stats import _update_history

    if not performance._has_dask():
        return _update_history(self._obj, "Optimization skipped (Dask not installed)")

    # Ensure data is lazy
    res = performance.apply_lazy_threshold(self._obj, threshold_mb=0.1)
    # Always calculate and apply recommended chunks for the target size
    recommendation = performance.get_chunk_recommendation(res, target_mb=target_mb)
    res = res.chunk(recommendation)

    return _update_history(res, f"Optimized for performance (target={target_mb}MB)")

peak_timing(dim='time')

Identify the coordinate value of the maximum.

Parameters

dim : str, optional Dimension along which to find peak. Default is 'time'.

Returns

xarray.Dataset Coordinate values.

Source code in src/monet_stats/accessor.py
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
def peak_timing(self, dim: str = "time") -> xr.Dataset:
    """
    Identify the coordinate value of the maximum.

    Parameters
    ----------
    dim : str, optional
        Dimension along which to find peak. Default is 'time'.

    Returns
    -------
    xarray.Dataset
        Coordinate values.
    """
    return analysis.peak_timing(self._obj, dim=dim)

percentile(q, dim='time', **kwargs)

Compute percentiles.

Parameters

q : float or list of float Percentile(s) (0-100). dim : str, optional Dimension over which to compute. Default is 'time'. **kwargs : Any Additional keyword arguments.

Returns

xarray.Dataset Computed percentiles.

Source code in src/monet_stats/accessor.py
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
def percentile(self, q: Union[float, List[float], np.ndarray], dim: str = "time", **kwargs: Any) -> xr.Dataset:
    """
    Compute percentiles.

    Parameters
    ----------
    q : float or list of float
        Percentile(s) (0-100).
    dim : str, optional
        Dimension over which to compute. Default is 'time'.
    **kwargs : Any
        Additional keyword arguments.

    Returns
    -------
    xarray.Dataset
        Computed percentiles.
    """
    return analysis.percentile(self._obj, q=q, dim=dim, **kwargs)

rechunk(chunks=None)

Apply new chunks to the Dataset (Aero Protocol provenance tracking).

Parameters

chunks : dict, optional New chunk sizes. If None, uses optimal recommendations (~100MB).

Returns

xarray.Dataset Rechunked Dataset.

Source code in src/monet_stats/accessor.py
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
def rechunk(self, chunks: Optional[dict] = None) -> xr.Dataset:
    """
    Apply new chunks to the Dataset (Aero Protocol provenance tracking).

    Parameters
    ----------
    chunks : dict, optional
        New chunk sizes. If None, uses optimal recommendations (~100MB).

    Returns
    -------
    xarray.Dataset
        Rechunked Dataset.
    """
    from .utils_stats import _update_history

    if not performance._has_dask():
        return _update_history(self._obj, "Rechunking skipped (Dask not installed)")

    if chunks is None:
        chunks = performance.get_chunk_recommendation(self._obj)

    res = self._obj.chunk(chunks)

    return _update_history(res, f"Rechunked with {chunks}")

resample_data(freq='MS', method='mean', dim='time', **kwargs)

Resample data to a new temporal frequency.

Parameters

freq : str, optional Resampling frequency. Default is 'MS'. method : str, optional Statistical method. Default is 'mean'. dim : str, optional Dimension along which to resample. Default is 'time'. **kwargs : Any Additional keyword arguments.

Returns

xarray.Dataset Resampled data.

Source code in src/monet_stats/accessor.py
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
def resample_data(self, freq: str = "MS", method: str = "mean", dim: str = "time", **kwargs: Any) -> xr.Dataset:
    """
    Resample data to a new temporal frequency.

    Parameters
    ----------
    freq : str, optional
        Resampling frequency. Default is 'MS'.
    method : str, optional
        Statistical method. Default is 'mean'.
    dim : str, optional
        Dimension along which to resample. Default is 'time'.
    **kwargs : Any
        Additional keyword arguments.

    Returns
    -------
    xarray.Dataset
        Resampled data.
    """
    return analysis.resample_data(self._obj, freq=freq, method=method, dim=dim, **kwargs)

rolling_mean_24h(dim='time', min_periods=18, center=True)

Compute rolling 24-hour mean.

Parameters

dim : str, optional Dimension along which to compute. Default is 'time'. min_periods : int, optional Minimum number of observations. Default is 18. center : bool, optional If True, center the labels. Default is True.

Returns

xarray.Dataset Rolling 24-hour mean.

Source code in src/monet_stats/accessor.py
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
def rolling_mean_24h(self, dim: str = "time", min_periods: int = 18, center: bool = True) -> xr.Dataset:
    """
    Compute rolling 24-hour mean.

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute. Default is 'time'.
    min_periods : int, optional
        Minimum number of observations. Default is 18.
    center : bool, optional
        If True, center the labels. Default is True.

    Returns
    -------
    xarray.Dataset
        Rolling 24-hour mean.
    """
    return analysis.rolling_mean_24h(self._obj, dim=dim, min_periods=min_periods, center=center)

rolling_mean_8h(dim='time', min_periods=6, center=True)

Compute rolling 8-hour mean.

Parameters

dim : str, optional Dimension along which to compute. Default is 'time'. min_periods : int, optional Minimum number of observations. Default is 6. center : bool, optional If True, center the labels. Default is True.

Returns

xarray.Dataset Rolling 8-hour mean.

Source code in src/monet_stats/accessor.py
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
def rolling_mean_8h(self, dim: str = "time", min_periods: int = 6, center: bool = True) -> xr.Dataset:
    """
    Compute rolling 8-hour mean.

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute. Default is 'time'.
    min_periods : int, optional
        Minimum number of observations. Default is 6.
    center : bool, optional
        If True, center the labels. Default is True.

    Returns
    -------
    xarray.Dataset
        Rolling 8-hour mean.
    """
    return analysis.rolling_mean_8h(self._obj, dim=dim, min_periods=min_periods, center=center)

seasonal_mean(dim='time', weighted=True)

Compute seasonal mean (DJF, MAM, JJA, SON).

Parameters

dim : str, optional Dimension along which to compute. Default is 'time'. weighted : bool, optional Weight by days in month. Default is True.

Returns

xarray.Dataset Seasonal means.

Source code in src/monet_stats/accessor.py
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
def seasonal_mean(self, dim: str = "time", weighted: bool = True) -> xr.Dataset:
    """
    Compute seasonal mean (DJF, MAM, JJA, SON).

    Parameters
    ----------
    dim : str, optional
        Dimension along which to compute. Default is 'time'.
    weighted : bool, optional
        Weight by days in month. Default is True.

    Returns
    -------
    xarray.Dataset
        Seasonal means.
    """
    return analysis.seasonal_mean(self._obj, dim=dim, weighted=weighted)

stats(obs_name='Obs', mod_name='Mod', threshold=0.0, minval=None, maxval=None)

Calculate summary statistics for observations and model results.

Parameters

obs_name : str, optional Name of observation variable. Default is 'Obs'. mod_name : str, optional Name of model variable. Default is 'Mod'. threshold : float, optional Threshold for contingency scores. Default is 0.0. minval : float, optional Minimum value for filtering. maxval : float, optional Maximum value for filtering.

Returns

dict Dictionary of calculated statistics.

Source code in src/monet_stats/accessor.py
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
def stats(
    self,
    obs_name: str = "Obs",
    mod_name: str = "Mod",
    threshold: float = 0.0,
    minval: Optional[float] = None,
    maxval: Optional[float] = None,
) -> dict:
    """
    Calculate summary statistics for observations and model results.

    Parameters
    ----------
    obs_name : str, optional
        Name of observation variable. Default is 'Obs'.
    mod_name : str, optional
        Name of model variable. Default is 'Mod'.
    threshold : float, optional
        Threshold for contingency scores. Default is 0.0.
    minval : float, optional
        Minimum value for filtering.
    maxval : float, optional
        Maximum value for filtering.

    Returns
    -------
    dict
        Dictionary of calculated statistics.
    """
    from . import stats

    return stats(
        self._obj,
        obs_name=obs_name,
        mod_name=mod_name,
        threshold=threshold,
        minval=minval,
        maxval=maxval,
    )

weighted_spatial_mean(lat_dim='lat', lon_dim='lon', weights=None)

Compute area-weighted spatial mean.

Parameters

lat_dim : str, optional Latitude dimension name. Default is 'lat'. lon_dim : str, optional Longitude dimension name. Default is 'lon'. weights : xarray.DataArray or numpy.ndarray, optional Custom weights.

Returns

xarray.Dataset Area-weighted spatial mean.

Source code in src/monet_stats/accessor.py
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
def weighted_spatial_mean(
    self,
    lat_dim: str = "lat",
    lon_dim: str = "lon",
    weights: Optional[Union[xr.DataArray, np.ndarray]] = None,
) -> xr.Dataset:
    """
    Compute area-weighted spatial mean.

    Parameters
    ----------
    lat_dim : str, optional
        Latitude dimension name. Default is 'lat'.
    lon_dim : str, optional
        Longitude dimension name. Default is 'lon'.
    weights : xarray.DataArray or numpy.ndarray, optional
        Custom weights.

    Returns
    -------
    xarray.Dataset
        Area-weighted spatial mean.
    """
    return analysis.weighted_spatial_mean(self._obj, lat_dim=lat_dim, lon_dim=lon_dim, weights=weights)