def count_level_crossings(low_data: xr.DataArray, high_data: xr.DataArray, threshold: str, freq: str) -> xr.DataArray: """Calculate the number of times low_data is below threshold while high_data is above threshold. First, the threshold is transformed to the same standard_name and units as the input data, then the thresholding is performed, and finally, the number of occurrences is counted. Parameters ---------- low_data: xr.DataArray Variable that must be under the threshold. high_data: xr.DataArray Variable that must be above the threshold. threshold: str Quantity. freq: str Resampling frequency. Returns ------- xarray.DataArray """ # Convert units to low_data high_data = convert_units_to(high_data, low_data) threshold = convert_units_to(threshold, low_data) lower = compare(low_data, "<", threshold) higher = compare(high_data, ">=", threshold) out = (lower & higher).resample(time=freq).sum() return to_agg_units(out, low_data, "count", dim="time")
def spell_length(data: xr.DataArray, threshold: str, condition: str, reducer: str, freq: str) -> xr.DataArray: """Calculate statistics on lengths of spells. First, the threshold is transformed to the same standard_name and units as the input data. Then the thresholding is performed as condition(data, threshold), i.e. if condition is <, data < threshold. Then the spells are determined, and finally the statistics according to the specified reducer are calculated. Parameters ---------- data : xr.DataArray threshold : str Quantity. condition : {">", "<", ">=", "<=", "==", "!="} Operator reducer : {'max', 'min', 'mean', 'sum'} Reducer. freq : str Resampling frequency. Returns ------- xarray.DataArray """ threshold = convert_units_to(threshold, data) cond = compare(data, condition, threshold) out = cond.resample(time=freq).map( rl.rle_statistics, reducer=reducer, dim="time", ) return to_agg_units(out, data, "count")
def count_occurrences(data: xr.DataArray, threshold: str, condition: str, freq: str) -> xr.DataArray: """Calculate the number of times some condition is met. First, the threshold is transformed to the same standard_name and units as the input data. Then the thresholding is performed as condition(data, threshold), i.e. if condition is `<`, then this counts the number of times `data < threshold`. Finally, count the number of occurrences when condition is met. Parameters ---------- data : xr.DataArray threshold : str Quantity. condition : {">", "<", ">=", "<=", "==", "!="} Operator. freq : str Resampling frequency. Returns ------- xarray.DataArray """ threshold = convert_units_to(threshold, data) cond = compare(data, condition, threshold) out = cond.resample(time=freq).sum() return to_agg_units(out, data, "count", dim="time")
def temperature_sum(data: xr.DataArray, threshold: str, condition: str, freq: str) -> xr.DataArray: """Calculate the temperature sum above/below a threshold. First, the threshold is transformed to the same standard_name and units as the input data. Then the thresholding is performed as condition(data, threshold), i.e. if condition is <, data < threshold. Finally, the sum is calculated for those data values that fulfill the condition after subtraction of the threshold value. If the sum is for values below the threshold the result is multiplied by -1. Parameters ---------- data : xr.DataArray threshold : str Quantity condition : {">", "<", ">=", "<=", "==", "!="} Operator freq : str Resampling frequency. Returns ------- xarray.DataArray """ threshold = convert_units_to(threshold, data) cond = compare(data, condition, threshold) direction = -1 if "<" in condition else 1 out = (data - threshold).where(cond).resample(time=freq).sum() out = direction * out return to_agg_units(out, data, "delta_prod")
def ice_days(tasmax: xarray.DataArray, thresh: str = "0 degC", freq: str = "YS"): # noqa: D401 r"""Number of ice/freezing days. Number of days where daily maximum temperatures are below 0℃. Parameters ---------- tasmax : xarray.DataArray Maximum daily temperature. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] Number of ice/freezing days. Notes ----- Let :math:`TX_{ij}` be the daily maximum temperature at day :math:`i` of period :math:`j`. Then counted is the number of days where: .. math:: TX_{ij} < 0℃ """ frz = convert_units_to(thresh, tasmax) out = threshold_count(tasmax, "<", frz, freq) return to_agg_units(out, tasmax, "count")
def frost_days(tasmin: xarray.DataArray, thresh: str = "0 degC", freq: str = "YS"): r"""Frost days index. Number of days where daily minimum temperatures are below 0℃. Parameters ---------- tasmin : xarray.DataArray Minimum daily temperature. thresh : str Freezing temperature. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] Frost days index. Notes ----- Let :math:`TN_{ij}` be the daily minimum temperature at day :math:`i` of period :math:`j`. Then counted is the number of days where: .. math:: TN_{ij} < 0℃ """ frz = convert_units_to(thresh, tasmin) out = threshold_count(tasmin, "<", frz, freq) return to_agg_units(out, tasmin, "count")
def degree_days(tas: xr.DataArray, thresh: str, condition: str) -> xr.DataArray: """Calculate the degree days below/above the temperature threshold. Parameters ---------- tas : xr.DataArray Mean daily temperature. thresh : str The temperature threshold. condition : {"<", ">"} Operator. Returns ------- xarray.DataArray """ thresh = convert_units_to(thresh, tas) if "<" in condition: out = (thresh - tas).clip(0) elif ">" in condition: out = (tas - thresh).clip(0) else: raise NotImplementedError(f"Condition not supported: '{condition}'.") out = to_agg_units(out, tas, op="delta_prod") return out
def rain_on_frozen_ground_days( pr: xarray.DataArray, tas: xarray.DataArray, thresh: str = "1 mm/d", freq: str = "YS", ) -> xarray.DataArray: # noqa: D401 """Number of rain on frozen ground events. Number of days with rain above a threshold after a series of seven days below freezing temperature. Precipitation is assumed to be rain when the temperature is above 0℃. Parameters ---------- pr : xarray.DataArray Mean daily precipitation flux. tas : xarray.DataArray Mean daily temperature. thresh : str Precipitation threshold to consider a day as a rain event. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] The number of rain on frozen ground events per period. Notes ----- Let :math:`PR_i` be the mean daily precipitation and :math:`TG_i` be the mean daily temperature of day :math:`i`. Then for a period :math:`j`, rain on frozen grounds days are counted where: .. math:: PR_{i} > Threshold [mm] and where .. math:: TG_{i} ≤ 0℃ is true for continuous periods where :math:`i ≥ 7` """ t = convert_units_to(thresh, pr) frz = convert_units_to("0 C", tas) def func(x, axis): """Check that temperature conditions are below 0 for seven days and above after.""" frozen = x == np.array([0, 0, 0, 0, 0, 0, 0, 1], bool) return frozen.all(axis=axis) tcond = (tas > frz).rolling(time=8).reduce(func) pcond = pr > t out = (tcond * pcond * 1).resample(time=freq).sum(dim="time") return to_agg_units(out, tas, "count")
def dry_spell_total_length( pr: xarray.DataArray, thresh: str = "1.0 mm", window: int = 3, op: str = "sum", freq: str = "YS", **indexer, ) -> xarray.DataArray: """ Total length of dry spells Total number of days in dry periods of a minimum length, during which the maximum or accumulated precipitation within a window of the same length is under a threshold. Parameters ---------- pr : xarray.DataArray Daily precipitation. thresh : str Accumulated precipitation value under which a period is considered dry. window : int Number of days where the maximum or accumulated precipitation is under threshold. op : {"max", "sum"} Reduce operation. freq : str Resampling frequency. indexer : Indexing parameters to compute the indicator on a temporal subset of the data. It accepts the same arguments as :py:func:`xclim.indices.generic.select_time`. Indexing is done after finding the dry days, but before finding the spells. Returns ------- xarray.DataArray The {freq} total number of days in dry periods of minimum {window} days. Notes ----- The algorithm assumes days before and after the timeseries are "wet", meaning that the condition for being considered part of a dry spell is stricter on the edges. For example, with `window=3` and `op='sum'`, the first day of the series is considered part of a dry spell only if the accumulated precipitation within the first 3 days is under the threshold. In comparison, a day in the middle of the series is considered part of a dry spell if any of the three 3-day periods of which it is part are considered dry (so a total of five days are included in the computation, compared to only 3.) """ pram = rate2amount(pr, out_units="mm") thresh = convert_units_to(thresh, pram) pram_pad = pram.pad(time=(0, window)) mask = getattr(pram_pad.rolling(time=window), op)() < thresh dry = (mask.rolling(time=window).sum() >= 1).shift(time=-(window - 1)) dry = dry.isel(time=slice(0, pram.time.size)).astype(float) out = select_time(dry, **indexer).resample(time=freq).sum("time") return to_agg_units(out, pram, "count")
def tx_tn_days_above( tasmin: xarray.DataArray, tasmax: xarray.DataArray, thresh_tasmin: str = "22 degC", thresh_tasmax: str = "30 degC", freq: str = "YS", ) -> xarray.DataArray: # noqa: D401 r"""Number of days with both hot maximum and minimum daily temperatures. The number of days per period with tasmin above a threshold and tasmax above another threshold. Parameters ---------- tasmin : xarray.DataArray Minimum daily temperature. tasmax : xarray.DataArray Maximum daily temperature. thresh_tasmin : str Threshold temperature for tasmin on which to base evaluation. thresh_tasmax : str Threshold temperature for tasmax on which to base evaluation. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] the number of days with tasmin > thresh_tasmin and tasmax > thresh_tasamax per period Notes ----- Let :math:`TX_{ij}` be the maximum temperature at day :math:`i` of period :math:`j`, :math:`TN_{ij}` the daily minimum temperature at day :math:`i` of period :math:`j`, :math:`TX_{thresh}` the threshold for maximum daily temperature, and :math:`TN_{thresh}` the threshold for minimum daily temperature. Then counted is the number of days where: .. math:: TX_{ij} > TX_{thresh} [℃] and where: .. math:: TN_{ij} > TN_{thresh} [℃] """ thresh_tasmax = convert_units_to(thresh_tasmax, tasmax) thresh_tasmin = convert_units_to(thresh_tasmin, tasmin) events = ((tasmin > thresh_tasmin) & (tasmax > thresh_tasmax)) * 1 out = events.resample(time=freq).sum(dim="time") return to_agg_units(out, tasmin, "count")
def warm_spell_duration_index( tasmax: xarray.DataArray, tx90: xarray.DataArray, window: int = 6, freq: str = "YS" ) -> xarray.DataArray: r"""Warm spell duration index. Number of days inside spells of a minimum number of consecutive days where the daily maximum temperature is above the 90th percentile. The 90th percentile should be computed for a 5-day moving window, centered on each calendar day in the 1961-1990 period. Parameters ---------- tasmax : xarray.DataArray Maximum daily temperature. tx90 : xarray.DataArray 90th percentile of daily maximum temperature. window : int Minimum number of days with temperature above threshold to qualify as a warm spell. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] Warm spell duration index. Examples -------- Note that this example does not use a proper 1961-1990 reference period. >>> from xclim.core.calendar import percentile_doy >>> from xclim.indices import warm_spell_duration_index >>> tasmax = xr.open_dataset(path_to_tasmax_file).tasmax.isel(lat=0, lon=0) >>> tx90 = percentile_doy(tasmax, per=90).sel(percentiles=90) >>> warm_spell_duration_index(tasmax, tx90) References ---------- From the Expert Team on Climate Change Detection, Monitoring and Indices (ETCCDMI). Used in Alexander, L. V., et al. (2006), Global observed changes in daily climate extremes of temperature and precipitation, J. Geophys. Res., 111, D05109, doi: 10.1029/2005JD006290. """ thresh = convert_units_to(tx90, tasmax) # Create time series out of doy values. thresh = resample_doy(thresh, tasmax) above = tasmax > thresh out = above.resample(time=freq).map( rl.windowed_run_count, window=window, dim="time" ) return to_agg_units(out, tasmax, "count")
def daily_freezethaw_cycles( tasmin: xarray.DataArray, tasmax: xarray.DataArray, thresh_tasmax: str = "0 degC", thresh_tasmin: str = "0 degC", freq: str = "YS", ) -> xarray.DataArray: # noqa: D401 r"""Number of days with a diurnal freeze-thaw cycle. The number of days where Tmax > thresh_tasmax and Tmin <= thresh_tasmin. Parameters ---------- tasmin : xarray.DataArray Minimum daily temperature. tasmax : xarray.DataArray Maximum daily temperature. thresh_tasmax : str The temperature threshold needed to trigger a thaw event. thresh_tasmin : str The temperature threshold needed to trigger a freeze event. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] Number of days with a diurnal freeze-thaw cycle Notes ----- Let :math:`TX_{i}` be the maximum temperature at day :math:`i` and :math:`TN_{i}` be the daily minimum temperature at day :math:`i`. Then the number of freeze thaw cycles during period :math:`\phi` is given by : .. math:: \sum_{i \in \phi} [ TX_{i} > 0℃ ] [ TN_{i} < 0℃ ] where :math:`[P]` is 1 if :math:`P` is true, and 0 if false. """ thaw_threshold = convert_units_to(thresh_tasmax, tasmax) freeze_threshold = convert_units_to(thresh_tasmin, tasmin) ft = (tasmin <= freeze_threshold) * (tasmax > thaw_threshold) * 1 out = ft.resample(time=freq).sum(dim="time") return to_agg_units(out, tasmin, "count")
def blowing_snow( snd: xarray.DataArray, sfcWind: xarray.DataArray, # noqa snd_thresh: str = "5 cm", sfcWind_thresh: str = "15 km/h", # noqa window: int = 3, freq: str = "AS-JUL", ) -> xarray.DataArray: """ Days with blowing snow events Number of days where both snowfall over the last days and daily wind speeds are above respective thresholds. Parameters ---------- snd : xarray.DataArray Surface snow depth. sfcWind : xr.DataArray Wind velocity snd_thresh : str Threshold on net snowfall accumulation over the last `window` days. sfcWind_thresh : str Wind speed threshold. window : int Period over which snow is accumulated before comparing against threshold. freq : str Resampling frequency. Returns ------- xarray.DataArray Number of days where snowfall and wind speeds are above respective thresholds. """ snd_thresh = convert_units_to(snd_thresh, snd) sfcWind_thresh = convert_units_to(sfcWind_thresh, sfcWind) # noqa # Net snow accumulation over the last `window` days snow = snd.diff(dim="time").rolling(time=window, center=False).sum() # Blowing snow conditions cond = (snow >= snd_thresh) * (sfcWind >= sfcWind_thresh) * 1 out = cond.resample(time=freq).sum(dim="time") out.attrs["units"] = to_agg_units(out, snd, "count") return out
def days_over_precip_thresh( pr: xarray.DataArray, per: xarray.DataArray, thresh: str = "1 mm/day", freq: str = "YS", ) -> xarray.DataArray: # noqa: D401 r"""Number of wet days with daily precipitation over a given percentile. Number of days over period where the precipitation is above a threshold defining wet days and above a given percentile for that day. Parameters ---------- pr : xarray.DataArray Mean daily precipitation flux. per : xarray.DataArray Daily percentile of wet day precipitation flux. thresh : str Precipitation value over which a day is considered wet. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] Count of days with daily precipitation above the given percentile [days]. Examples -------- >>> from xclim.indices import days_over_precip_thresh >>> pr = xr.open_dataset(path_to_pr_file).pr >>> p75 = pr.quantile(.75, dim="time", keep_attrs=True) >>> r75p = days_over_precip_thresh(pr, p75) """ per = convert_units_to(per, pr) thresh = convert_units_to(thresh, pr) tp = np.maximum(per, thresh) if "dayofyear" in per.coords: # Create time series out of doy values. tp = resample_doy(tp, pr) # Compute the days where precip is both over the wet day threshold and the percentile threshold. out = threshold_count(pr, ">", tp, freq) return to_agg_units(out, pr, "count")
def heat_wave_total_length( tasmin: xarray.DataArray, tasmax: xarray.DataArray, thresh_tasmin: str = "22.0 degC", thresh_tasmax: str = "30 degC", window: int = 3, freq: str = "YS", ) -> xarray.DataArray: r"""Heat wave total length. Total length of heat waves over a given period. A heat wave is defined as an event where the minimum and maximum daily temperature both exceeds specific thresholds over a minimum number of days. This the sum of all days in such events. Parameters ---------- tasmin : xarray.DataArray Minimum daily temperature. tasmax : xarray.DataArray Maximum daily temperature. thresh_tasmin : str The minimum temperature threshold needed to trigger a heatwave event. thresh_tasmax : str The maximum temperature threshold needed to trigger a heatwave event. window : int Minimum number of days with temperatures above thresholds to qualify as a heatwave. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] Total length of heatwave at the requested frequency. Notes ----- See notes and references of `heat_wave_max_length` """ thresh_tasmax = convert_units_to(thresh_tasmax, tasmax) thresh_tasmin = convert_units_to(thresh_tasmin, tasmin) cond = (tasmin > thresh_tasmin) & (tasmax > thresh_tasmax) out = cond.resample(time=freq).map(rl.windowed_run_count, window=window) return to_agg_units(out, tasmin, "count")
def high_precip_low_temp( pr: xarray.DataArray, tas: xarray.DataArray, pr_thresh: str = "0.4 mm/d", tas_thresh: str = "-0.2 degC", freq: str = "YS", ) -> xarray.DataArray: # noqa: D401 """Number of days with precipitation above threshold and temperature below threshold. Number of days where precipitation is greater or equal to some threshold, and temperatures are colder than some threshold. This can be used for example to identify days with the potential for freezing rain or icing conditions. Parameters ---------- pr : xarray.DataArray Mean daily precipitation flux. tas : xarray.DataArray Daily mean, minimum or maximum temperature. pr_thresh : str Precipitation threshold to exceed. tas_thresh : str Temperature threshold not to exceed. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] Count of days with high precipitation and low temperatures. Example ------- To compute the number of days with intense rainfall while minimum temperatures dip below -0.2C: >>> pr = xr.open_dataset(path_to_pr_file).pr >>> tasmin = xr.open_dataset(path_to_tasmin_file).tasmin >>> high_precip_low_temp(pr, tas=tasmin, pr_thresh="10 mm/d", tas_thresh="-0.2 degC") """ pr_thresh = convert_units_to(pr_thresh, pr) tas_thresh = convert_units_to(tas_thresh, tas) cond = (pr >= pr_thresh) * (tas < tas_thresh) * 1 out = cond.resample(time=freq).sum(dim="time") return to_agg_units(out, pr, "count")
def tn90p( tasmin: xarray.DataArray, t90: xarray.DataArray, freq: str = "YS" ) -> xarray.DataArray: # noqa: D401 r"""Number of days with daily minimum temperature over the 90th percentile. Number of days with daily minimum temperature over the 90th percentile. Parameters ---------- tasmin : xarray.DataArray Minimum daily temperature. t90 : xarray.DataArray 90th percentile of daily minimum temperature. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] Count of days with daily minimum temperature below the 10th percentile [days]. Notes ----- The 90th percentile should be computed for a 5 day window centered on each calendar day for a reference period. Examples -------- >>> from xclim.core.calendar import percentile_doy >>> from xclim.indices import tn90p >>> tas = xr.open_dataset(path_to_tas_file).tas >>> t90 = percentile_doy(tas, per=90).sel(percentiles=90) >>> hot_days = tn90p(tas, t90) """ t90 = convert_units_to(t90, tasmin) # Create time series out of doy values. thresh = resample_doy(t90, tasmin) # Identify the days with min temp above 90th percentile. out = threshold_count(tasmin, ">", thresh, freq) return to_agg_units(out, tasmin, "count")
def cold_spell_duration_index( tasmin: xarray.DataArray, tn10: xarray.DataArray, window: int = 6, freq: str = "YS" ) -> xarray.DataArray: r"""Cold spell duration index. Number of days with at least six consecutive days where the daily minimum temperature is below the 10th percentile. Parameters ---------- tasmin : xarray.DataArray Minimum daily temperature. tn10 : xarray.DataArray 10th percentile of daily minimum temperature with `dayofyear` coordinate. window : int Minimum number of days with temperature below threshold to qualify as a cold spell. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] Count of days with at least six consecutive days where the daily minimum temperature is below the 10th percentile. Notes ----- Let :math:`TN_i` be the minimum daily temperature for the day of the year :math:`i` and :math:`TN10_i` the 10th percentile of the minimum daily temperature over the 1961-1990 period for day of the year :math:`i`, the cold spell duration index over period :math:`\phi` is defined as: .. math:: \sum_{i \in \phi} \prod_{j=i}^{i+6} \left[ TN_j < TN10_j \right] where :math:`[P]` is 1 if :math:`P` is true, and 0 if false. References ---------- From the Expert Team on Climate Change Detection, Monitoring and Indices (ETCCDMI). Examples -------- # Note that this example does not use a proper 1961-1990 reference period. >>> from xclim.core.calendar import percentile_doy >>> from xclim.indices import cold_spell_duration_index >>> tasmin = xr.open_dataset(path_to_tasmin_file).tasmin.isel(lat=0, lon=0) >>> tn10 = percentile_doy(tasmin, per=10).sel(percentiles=10) >>> cold_spell_duration_index(tasmin, tn10) """ tn10 = convert_units_to(tn10, tasmin) # Create time series out of doy values. thresh = resample_doy(tn10, tasmin) below = tasmin < thresh out = below.resample(time=freq).map( rl.windowed_run_count, window=window, dim="time" ) return to_agg_units(out, tasmin, "count")
def heat_wave_max_length( tasmin: xarray.DataArray, tasmax: xarray.DataArray, thresh_tasmin: str = "22.0 degC", thresh_tasmax: str = "30 degC", window: int = 3, freq: str = "YS", ) -> xarray.DataArray: r"""Heat wave max length. Maximum length of heat waves over a given period. A heat wave is defined as an event where the minimum and maximum daily temperature both exceeds specific thresholds over a minimum number of days. By definition heat_wave_max_length must be >= window. Parameters ---------- tasmin : xarray.DataArray Minimum daily temperature. tasmax : xarray.DataArray Maximum daily temperature. thresh_tasmin : str The minimum temperature threshold needed to trigger a heatwave event. thresh_tasmax : str The maximum temperature threshold needed to trigger a heatwave event. window : int Minimum number of days with temperatures above thresholds to qualify as a heatwave. freq : str Resampling frequency. Returns ------- xarray.DataArray, [time] Maximum length of heatwave at the requested frequency. Notes ----- The thresholds of 22° and 25°C for night temperatures and 30° and 35°C for day temperatures were selected by Health Canada professionals, following a temperature–mortality analysis. These absolute temperature thresholds characterize the occurrence of hot weather events that can result in adverse health outcomes for Canadian communities (Casati et al., 2013). In Robinson (2001), the parameters would be `thresh_tasmin=27.22, thresh_tasmax=39.44, window=2` (81F, 103F). References ---------- Casati, B., A. Yagouti, and D. Chaumont, 2013: Regional Climate Projections of Extreme Heat Events in Nine Pilot Canadian Communities for Public Health Planning. J. Appl. Meteor. Climatol., 52, 2669–2698, https://doi.org/10.1175/JAMC-D-12-0341.1 Robinson, P.J., 2001: On the Definition of a Heat Wave. J. Appl. Meteor., 40, 762–775, https://doi.org/10.1175/1520-0450(2001)040<0762:OTDOAH>2.0.CO;2 """ thresh_tasmax = convert_units_to(thresh_tasmax, tasmax) thresh_tasmin = convert_units_to(thresh_tasmin, tasmin) cond = (tasmin > thresh_tasmin) & (tasmax > thresh_tasmax) max_l = cond.resample(time=freq).map(rl.longest_run, dim="time") out = max_l.where(max_l >= window, 0) return to_agg_units(out, tasmax, "count")
def effective_growing_degree_days( tasmax: xarray.DataArray, tasmin: xarray.DataArray, *, thresh: str = "5 degC", method: str = "bootsma", after_date: DayOfYearStr = "07-01", dim: str = "time", freq: str = "YS", ) -> xarray.DataArray: r"""Effective growing degree days. Growing degree days based on a dynamic start and end of the growing season. Parameters ---------- tasmax: xr.DataArray Daily mean temperature. tasmin: xr.DataArray Daily minimum temperature. thresh: str The minimum temperature threshold. method: {"bootsma", "qian"} The window method used to determine the temperature-based start date. For "bootsma", the start date is defined as 10 days after the average temperature exceeds a threshold (5 degC). For "qian", the start date is based on a weighted 5-day rolling average, based on `qian_weighted_mean_average()`. after_date : str Date of the year after which to look for the first frost event. Should have the format '%m-%d'. dim: str Time dimension. freq : str Resampling frequency. Returns ------- xarray.DataArray Notes ----- The effective growing degree days for a given year :math:`EGDD_i` can be calculated as follows: .. math:: EGDD_i = \sum_{i=\text{j_{start}}^{\text{j_{end}}} max\left(TG - Thresh, 0 \right) Where :math:`TG` is the mean daly temperature, and :math:`j_{start}` and :math:`j_{end}` are the start and end dates of the growing season. The growing season start date methodology is determined via the `method` flag. For "bootsma", the start date is defined as 10 days after the average temperature exceeds a threshold (5 degC). For "qian", the start date is based on a weighted 5-day rolling average, based on `qian_weighted_mean_average()`. The end date is determined as the day preceding the first day with minimum temperature below 0 degC. References ---------- Indice originally published in Bootsma, A., & Gameda and D.W. McKenney, S. (2005). Impacts of potential climate change on selected agroclimatic indices in Atlantic Canada. Canadian Journal of Soil Science, 85(2), 329‑343. https://doi.org/10.4141/S04-019 """ tasmax = convert_units_to(tasmax, "degC") tasmin = convert_units_to(tasmin, "degC") thresh = convert_units_to(thresh, "degC") tas = (tasmin + tasmax) / 2 tas.attrs["units"] = "degC" if method.lower() == "bootsma": fda = first_day_above(tasmin=tas, thresh="5.0 degC", window=1, freq=freq) start = fda + 10 elif method.lower() == "qian": tas_weighted = qian_weighted_mean_average(tas=tas, dim=dim) start = freshet_start(tas_weighted, thresh=thresh, window=5, freq=freq) else: raise NotImplementedError(f"Method: {method}.") # The day before the first day below 0 degC end = (first_day_below( tasmin=tasmin, thresh="0 degC", after_date=after_date, window=1, freq=freq, ) - 1) deg_days = (tas - thresh).clip(min=0) egdd = aggregate_between_dates(deg_days, start=start, end=end, freq=freq) return to_agg_units(egdd, tas, op="delta_prod")