def test_datetime_to_decimal_year(source_cal, exp180): times = xr.DataArray( date_range( "2004-01-01", "2004-12-30", freq="D", calendar=source_cal or "default" ), dims=("time",), name="time", ) decy = datetime_to_decimal_year(times, calendar=source_cal) np.testing.assert_almost_equal(decy[180] - 2004, exp180)
def potential_evapotranspiration( tasmin: Optional[xr.DataArray] = None, tasmax: Optional[xr.DataArray] = None, tas: Optional[xr.DataArray] = None, method: str = "BR65", peta: Optional[float] = 0.00516409319477, petb: Optional[float] = 0.0874972822289, ) -> xr.DataArray: """Potential evapotranspiration. The potential for water evaporation from soil and transpiration by plants if the water supply is sufficient, according to a given method. Parameters ---------- tasmin : xarray.DataArray Minimum daily temperature. tasmax : xarray.DataArray Maximum daily temperature. tas : xarray.DataArray Mean daily temperature. method : {"baierrobertson65", "BR65", "hargreaves85", "HG85", "thornthwaite48", "TW48", "mcguinnessbordne05", "MB05"} Which method to use, see notes. peta : float Used only with method MB05 as :math:`a` for calculation of PET, see Notes section. Default value resulted from calibration of PET over the UK. petb : float Used only with method MB05 as :math:`b` for calculation of PET, see Notes section. Default value resulted from calibration of PET over the UK. Returns ------- xarray.DataArray Notes ----- Available methods are: - "baierrobertson65" or "BR65", based on [baierrobertson65]_. Requires tasmin and tasmax, daily [D] freq. - "hargreaves85" or "HG85", based on [hargreaves85]_. Requires tasmin and tasmax, daily [D] freq. (optional: tas can be given in addition of tasmin and tasmax). - "mcguinnessbordne05" or "MB05", based on [tanguy2018]_. Requires tas, daily [D] freq, with latitudes 'lat'. - "thornthwaite48" or "TW48", based on [thornthwaite48]_. Requires tasmin and tasmax, monthly [MS] or daily [D] freq. (optional: tas can be given instead of tasmin and tasmax). The McGuinness-Bordne [McGuinness1972]_ equation is: .. math:: PET[mm day^{-1}] = a * \frac{S_0}{\\lambda}T_a + b *\frsc{S_0}{\\lambda} where :math:`a` and :math:`b` are empirical parameters; :math:`S_0` is the extraterrestrial radiation [MJ m-2 day-1]; :math:`\\lambda` is the latent heat of vaporisation [MJ kg-1] and :math:`T_a` is the air temperature [°C]. The equation was originally derived for the USA, with :math:`a=0.0147` and :math:`b=0.07353`. The default parameters used here are calibrated for the UK, using the method described in [Tanguy2018]_. References ---------- .. [baierrobertson65] Baier, W., & Robertson, G. W. (1965). Estimation of latent evaporation from simple weather observations. Canadian journal of plant science, 45(3), 276-284. .. [hargreaves85] Hargreaves, G. H., & Samani, Z. A. (1985). Reference crop evapotranspiration from temperature. Applied engineering in agriculture, 1(2), 96-99. .. [tanguy2018] Tanguy, M., Prudhomme, C., Smith, K., & Hannaford, J. (2018). Historical gridded reconstruction of potential evapotranspiration for the UK. Earth System Science Data, 10(2), 951-968. .. [McGuinness1972] McGuinness, J. L., & Bordne, E. F. (1972). A comparison of lysimeter-derived potential evapotranspiration with computed values (No. 1452). US Department of Agriculture. .. [thornthwaite48] Thornthwaite, C. W. (1948). An approach toward a rational classification of climate. Geographical review, 38(1), 55-94. """ if method in ["baierrobertson65", "BR65"]: tasmin = convert_units_to(tasmin, "degF") tasmax = convert_units_to(tasmax, "degF") latr = (tasmin.lat * np.pi) / 180 gsc = 0.082 # MJ/m2/min # julian day fraction jd_frac = (datetime_to_decimal_year(tasmin.time) % 1) * 2 * np.pi ds = 0.409 * np.sin(jd_frac - 1.39) dr = 1 + 0.033 * np.cos(jd_frac) omega = np.arccos(-np.tan(latr) * np.tan(ds)) re = ((24 * 60 / np.pi) * gsc * dr * (omega * np.sin(latr) * np.sin(ds) + np.cos(latr) * np.cos(ds) * np.sin(omega))) # MJ/m2/day re = re / 4.1864e-2 # cal/cm2/day # Baier et Robertson(1965) formula out = 0.094 * (-87.03 + 0.928 * tasmax + 0.933 * (tasmax - tasmin) + 0.0486 * re) out = out.clip(0) elif method in ["hargreaves85", "HG85"]: tasmin = convert_units_to(tasmin, "degC") tasmax = convert_units_to(tasmax, "degC") if tas is None: tas = (tasmin + tasmax) / 2 else: tas = convert_units_to(tas, "degC") latr = (tasmin.lat * np.pi) / 180 gsc = 0.082 # MJ/m2/min lv = 2.5 # MJ/kg # julian day fraction jd_frac = (datetime_to_decimal_year(tasmin.time) % 1) * 2 * np.pi ds = 0.409 * np.sin(jd_frac - 1.39) dr = 1 + 0.033 * np.cos(jd_frac) omega = np.arccos(-np.tan(latr) * np.tan(ds)) ra = ((24 * 60 / np.pi) * gsc * dr * (omega * np.sin(latr) * np.sin(ds) + np.cos(latr) * np.cos(ds) * np.sin(omega))) # MJ/m2/day # Hargreaves and Samani(1985) formula out = (0.0023 * ra * (tas + 17.8) * (tasmax - tasmin)**0.5) / lv out = out.clip(0) elif method in ["mcguinnessbordne05", "MB05"]: if tas is None: tasmin = convert_units_to(tasmin, "degC") tasmax = convert_units_to(tasmax, "degC") tas = (tasmin + tasmax) / 2 tas.attrs["units"] = "degC" tas = convert_units_to(tas, "degC") tasK = convert_units_to(tas, "K") latr = (tas.lat * np.pi) / 180 jd_frac = (datetime_to_decimal_year(tas.time) % 1) * 2 * np.pi S = 1367.0 # Set solar constant [W/m2] ds = 0.409 * np.sin(jd_frac - 1.39) # solar declination ds [radians] omega = np.arccos(-np.tan(latr) * np.tan(ds)) # sunset hour angle [radians] dr = 1.0 + 0.03344 * np.cos( jd_frac - 0.048869) # Calculate relative distance to sun ext_rad = (S * 86400 / np.pi * dr * (omega * np.sin(ds) * np.sin(latr) + np.sin(omega) * np.cos(ds) * np.cos(latr))) latentH = 4185.5 * (751.78 - 0.5655 * tasK) radDIVlat = ext_rad / latentH # parameters from calibration provided by Dr Maliko Tanguy @ CEH # (calibrated for PET over the UK) a = peta b = petb out = radDIVlat * a * tas + radDIVlat * b elif method in ["thornthwaite48", "TW48"]: if tas is None: tasmin = convert_units_to(tasmin, "degC") tasmax = convert_units_to(tasmax, "degC") tas = (tasmin + tasmax) / 2 else: tas = convert_units_to(tas, "degC") tas = tas.clip(0) tas = tas.resample(time="MS").mean(dim="time") latr = (tas.lat * np.pi) / 180 # rad start = "-".join([ str(tas.time[0].dt.year.values), f"{tas.time[0].dt.month.values:02d}", "01", ]) end = "-".join([ str(tas.time[-1].dt.year.values), f"{tas.time[-1].dt.month.values:02d}", str(tas.time[-1].dt.daysinmonth.values), ]) time_v = xr.DataArray( date_range(start, end, freq="D", calendar="standard"), dims="time", name="time", ) # julian day fraction jd_frac = (datetime_to_decimal_year(time_v) % 1) * 2 * np.pi ds = 0.409 * np.sin(jd_frac - 1.39) omega = np.arccos(-np.tan(latr) * np.tan(ds)) * 180 / np.pi # degrees # monthly-mean daytime length (multiples of 12 hours) dl = 2 * omega / (15 * 12) dl_m = dl.resample(time="MS").mean(dim="time") # annual heat index id_m = (tas / 5)**1.514 id_y = id_m.resample(time="YS").sum(dim="time") tas_idy_a = [] for base_time, indexes in tas.resample(time="YS").groups.items(): tas_y = tas.isel(time=indexes) id_v = id_y.sel(time=base_time) a = 6.75e-7 * id_v**3 - 7.71e-5 * id_v**2 + 0.01791 * id_v + 0.49239 frac = (10 * tas_y / id_v)**a tas_idy_a.append(frac) tas_idy_a = xr.concat(tas_idy_a, dim="time") # Thornthwaite(1948) formula out = 1.6 * dl_m * tas_idy_a # cm/month out = 10 * out # mm/month else: raise NotImplementedError(f"'{method}' method is not implemented.") out.attrs["units"] = "mm" return amount2rate(out, out_units="kg m-2 s-1")