pydaymet.pydaymet
Contents
pydaymet.pydaymet
#
Access the Daymet database for both single single pixel and gridded queries.
Module Contents#
- pydaymet.pydaymet.get_bycoords(coords, dates, coords_id=None, crs=4326, variables=None, region='na', time_scale='daily', pet=None, pet_params=None, snow=False, snow_params=None, ssl=True, to_xarray=False)#
Get point-data from the Daymet database at 1-km resolution.
This function uses THREDDS data service to get the coordinates and supports getting monthly and annual summaries of the climate data directly from the server.
- Parameters
coords (
tuple
orlist
oftuples
) – Coordinates of the location(s) of interest as a tuple (x, y)dates (
tuple
orlist
, optional) – Start and end dates as a tuple (start, end) or a list of years[2001, 2010, ...]
.coords_id (
list
ofint
orstr
, optional) – A list of identifiers for the coordinates. This option only applies whento_xarray
is set toTrue
. If not provided, the coordinates will be enumerated.crs (
str
,int
, orpyproj.CRS
, optional) – The CRS of the input coordinates, defaults toEPSG:4326
.variables (
str
orlist
) – List of variables to be downloaded. The acceptable variables are:tmin
,tmax
,prcp
,srad
,vp
,swe
,dayl
Descriptions can be found here.region (
str
, optional) – Target region in the US, defaults tona
. Acceptable values are:na
: Continental North Americahi
: Hawaiipr
: Puerto Rico
time_scale (
str
, optional) – Data time scale which can bedaily
,monthly
(monthly summaries), orannual
(annual summaries). Defaults todaily
.pet (
str
, optional) – Method for computing PET. Supported methods arepenman_monteith
,priestley_taylor
,hargreaves_samani
, and None (don’t compute PET). Thepenman_monteith
method is based on Allen et al.1 assuming that soil heat flux density is zero. Thepriestley_taylor
method is based on Priestley and TAYLOR2 assuming that soil heat flux density is zero. Thehargreaves_samani
method is based on Hargreaves and Samani3. Defaults toNone
.pet_params (
dict
, optional) – Model-specific parameters as a dictionary that is passed to the PET function. Defaults toNone
.snow (
bool
, optional) – Compute snowfall from precipitation and minimum temperature. Defaults toFalse
.snow_params (
dict
, optional) – Model-specific parameters as a dictionary that is passed to the snowfall function. These parameters are only used ifsnow
isTrue
. Two parameters are required:t_rain
(deg C) which is the threshold for temperature for considering rain andt_snow
(deg C) which is the threshold for temperature for considering snow. The default values are{'t_rain': 2.5, 't_snow': 0.6}
that are adopted from https://doi.org/10.5194/gmd-11-1077-2018.ssl (
bool
, optional) – Whether to verify SSL certification, defaults toTrue
.to_xarray (
bool
, optional) – Return the data as anxarray.Dataset
. Defaults toFalse
.
- Returns
pandas.DataFrame
orxarray.Dataset
– Daily climate data for a single or list of locations.
Examples
>>> import pydaymet as daymet >>> coords = (-1431147.7928, 318483.4618) >>> dates = ("2000-01-01", "2000-12-31") >>> clm = daymet.get_bycoords( ... coords, ... dates, ... crs="epsg:3542", ... pet="hargreaves_samani", ... ) >>> clm["pet (mm/day)"].mean() 3.713
References
- 1(1,2)
Richard G Allen, Luis S Pereira, Dirk Raes, Martin Smith, and others. Crop evapotranspiration-guidelines for computing crop water requirements-fao irrigation and drainage paper 56. Fao, Rome, 300(9):D05109, 1998.
- 2(1,2)
Charles Henry Brian Priestley and Robert Joseph TAYLOR. On the assessment of surface heat flux and evaporation using large-scale parameters. Monthly weather review, 100(2):81–92, 1972.
- 3(1,2)
George H. Hargreaves and Zohrab A. Samani. Estimating potential evapotranspiration. Journal of the Irrigation and Drainage Division, 108(3):225–230, sep 1982. URL: https://doi.org/10.1061%2Fjrcea4.0001390, doi:10.1061/jrcea4.0001390.
- pydaymet.pydaymet.get_bygeom(geometry, dates, crs=4326, variables=None, region='na', time_scale='daily', pet=None, pet_params=None, snow=False, snow_params=None, ssl=True)#
Get gridded data from the Daymet database at 1-km resolution.
- Parameters
geometry (
Polygon
,MultiPolygon
, orbbox
) – The geometry of the region of interest.dates (
tuple
orlist
, optional) – Start and end dates as a tuple (start, end) or a list of years [2001, 2010, …].crs (
str
,int
, orpyproj.CRS
, optional) – The CRS of the input geometry, defaults to epsg:4326.variables (
str
orlist
) – List of variables to be downloaded. The acceptable variables are:tmin
,tmax
,prcp
,srad
,vp
,swe
,dayl
Descriptions can be found here.region (
str
, optional) – Region in the US, defaults to na. Acceptable values are:na: Continental North America
hi: Hawaii
pr: Puerto Rico
time_scale (
str
, optional) – Data time scale which can be daily, monthly (monthly average), or annual (annual average). Defaults to daily.pet (
str
, optional) – Method for computing PET. Supported methods arepenman_monteith
,priestley_taylor
,hargreaves_samani
, and None (don’t compute PET). Thepenman_monteith
method is based on Allen et al.1 assuming that soil heat flux density is zero. Thepriestley_taylor
method is based on Priestley and TAYLOR2 assuming that soil heat flux density is zero. Thehargreaves_samani
method is based on Hargreaves and Samani3. Defaults toNone
.pet_params (
dict
, optional) – Model-specific parameters as a dictionary that is passed to the PET function. Defaults toNone
.snow (
bool
, optional) – Compute snowfall from precipitation and minimum temperature. Defaults toFalse
.snow_params (
dict
, optional) – Model-specific parameters as a dictionary that is passed to the snowfall function. These parameters are only used ifsnow
isTrue
. Two parameters are required:t_rain
(deg C) which is the threshold for temperature for considering rain andt_snow
(deg C) which is the threshold for temperature for considering snow. The default values are{'t_rain': 2.5, 't_snow': 0.6}
that are adopted from https://doi.org/10.5194/gmd-11-1077-2018.ssl (
bool
, optional) – Whether to verify SSL certification, defaults toTrue
.
- Returns
xarray.Dataset
– Daily climate data within the target geometry.
Examples
>>> from shapely.geometry import Polygon >>> import pydaymet as daymet >>> geometry = Polygon( ... [[-69.77, 45.07], [-69.31, 45.07], [-69.31, 45.45], [-69.77, 45.45], [-69.77, 45.07]] ... ) >>> clm = daymet.get_bygeom(geometry, 2010, variables="tmin", time_scale="annual") >>> clm["tmin"].mean().compute().item() 1.361
References