View source: R/extractCatchmentData.R
extractCatchmentData | R Documentation |
extractCatchmentData
extracts catchment average climate data from netCDF files containing Australian climate data.extractCatchmentData extracts the AWAP climate data for each point or polygon. For the latter, either the daily spatial mean and variance (or user defined function) of each climate metric is calculated or the spatial data is returned.
extractCatchmentData(
ncdfFilename = file.path(getwd(), "AWAP.nc"),
ncdfSolarFilename = file.path(getwd(), "AWAP_solar.nc"),
extractFrom = as.Date("1900-01-01", "%Y-%m-%d"),
extractTo = as.Date(Sys.Date(), "%Y-%m-%d"),
getPrecip = TRUE,
getTmin = TRUE,
getTmax = TRUE,
getVprp = TRUE,
getSolarrad = TRUE,
getET = TRUE,
DEM = "",
catchments = "",
temporal.timestep = "daily",
temporal.function.name = "mean",
spatial.function.name = "var",
interpMethod = "",
ET.function = "ET.MortonCRAE",
ET.Mortons.est = "potential ET",
ET.Turc.humid = F,
ET.timestep = "monthly",
ET.interp_missing_days = T,
ET.interp_missing_entries = T,
ET.interp_abnormal = T,
ET.constants = list()
)
ncdfFilename |
is a full file name (as string) to the netCDF file. |
ncdfSolarFilename |
is the full file name (as string) to the netCDF file. |
extractFrom |
is a date string specifying the start date for data extraction. The default is |
extractTo |
is a date string specifying the end date for the data extraction. The default is today's date as YYYY-MM-DD. |
getPrecip |
logical variable for extracting precipitation. Default is |
getTmin |
logical variable for extracting Tmin. Default is |
getTmax |
logical variable for extracting Tmax. Default is |
getVprp |
logical variable for extracting vapour pressure. Default is |
getSolarrad |
logical variable for extracting solar radiation. Default is |
getET |
logical variable for calculating Morton's potential ET. Note, to calculate set |
DEM |
is either the full file name to a ESRI ASCII grid (as lat/long and using GDA94) or a raster class grid object. The DEM is used
for the calculation of Morton's PET. The Australian 9 second DEM can be loaded using |
catchments |
is either the full file name to an ESRI shape file of points or polygons (latter assumed to be catchment boundaries) or a shape file already imported using readShapeSpatial(). Either way the shape file must be in long/lat (i.e. not projected), use the ellipsoid GRS 80, and the first column should be a unique ID. |
temporal.timestep |
character string for the time step of the output data. The options are |
temporal.function.name |
character string for the function name applied to aggregate the daily data to |
spatial.function.name |
character string for the function name applied to estimate the daily spatial spread in each variable. If |
interpMethod |
character string defining the method for interpolating the gridded data (see |
ET.function |
character string for the evapotranspiration function to be used. The methods that can be derived from the AWAP data are are |
ET.Mortons.est |
character string for the type of Mortons Et estimate. For |
ET.Turc.humid |
logical variable for the Turc function using the humid adjustment.See |
ET.timestep |
character string for the evapotranpiration time step. Options are |
ET.interp_missing_days |
T or F, indicating if missing days should be interpolated for PET calculation. Default is |
ET.interp_missing_entries |
T or F, indicating if missing data entries should be interpolated for PET calculation. Default is |
ET.interp_abnormal |
T or F, indicating if abnormal valuses should be interpolated for PET calculation. Default is |
ET.constants |
list of constants from Evapotranspiration package required for ET calculations. To get the data use the command |
Daily data is extracted and can be aggregated to a weekly, monthly, quarterly, annual or a user-defined timestep using a user-defined funcion
(e.g. sum, mean, min, max as defined by temporal.function.name
). The temporally aggreated data at each grid cell is then used to derive the spatial
mean or the spatial variance (or any other function as defined by spatial.function.name
).
The calculation of the spatial mean uses the fraction of each AWAP grid cell within the catchment polygon. The variance calculation (or user defined function) does not use the fraction of the grid cell and returns NA if there are <2 grid cells in the catchment boundary. Prior to the spatial aggregation, evapotranspiration (ET) can also calculated; after which, say, the mean and variance PET can be calculated.
The data extraction will by default be undertaken from 1/1/1900 to yesterday, even if the netCDF grids were only built for a subset of this time period. If the latter situation applies, it is recommended that the extraction start and end dates are input by the user.
The ET can be calculated using one of eight methods at a user defined calculation time-step; that is the ET.timestep
defines the
time step at which the estimates are serived and differs from the output timestep as defined by temporal.function.name
). When ET.timestep
is monthly or annual then
the ET estimate is linearly interpolated to a daily time step (using zoo:na.spline()) and then constrained to >=0. In calculating ET, the input data
is pre-processed using Evapotranspiration::ReadInputs() such that missing days, missing enteries and abnormal values are interpolated
(by default) with the former two interpolated using the "DoY average", i.e. replacement with same day-of-the-year average. Additionally, when AWAP solar
radiation is required for the ET function, data is only available from 1/1/1990. To derive ET values <1990, the average solar radiation for each day of the year from
1/1/990 to "extractTo" is derived (i.e. 365 values) and then applied to each day prior to 1990. Importantly, in this situation the estimates of ET <1990
are dependent upon the end date extracted. Re-running the estimation of ET with a later extractTo data will change the estimates of ET
prior to 1990.
Also, when "catchments" is points (not polygons), then the netCDF grids are interpolate using bilinear interpolation of the closest 4 grid cells.
Lastly, data is extracted for all time points and no temporal infilling is undertaken if the grid cells are blank.
When catchments
are polygons and spatial.function.name
is not NA
or ""
, then the returned variable is a list variable containing two data.frames. The first is the areal aggreggated climate
metrics named catchmentTemporal.
with a suffix as defined by temporal.function.name
). The second is the measure of spatial variability
named catchmentSpatial.
with a suffix as defined by spatial.function.name
).
When catchments
are polygons and spatial.function.name
does equal NA
or ""
, then the returned variable is a sp::SpatialPixelsDataFrame
where the first colum is the catchment IDs
and the latter columns are the results for each variable at each time point as defined by temporal.timestep
.
When catchments
are points, the returned variable is a data.frame containing daily climate data at each point.
makeNetCDF_file
for building the NetCDF files of daily climate data.
# The example shows how to extract and save data.
# For an additional example see \url{https://github.com/peterson-tim-j/AWAPer/blob/master/README.md}
#---------------------------------------
library(sp)
# Set dates for building netCDFs and extracting data.
# Note, to reduce runtime this is done only a fortnight (14 days).
startDate = as.Date("2000-01-01","%Y-%m-%d")
endDate = as.Date("2000-01-14","%Y-%m-%d")
# Set names for netCDF files.
ncdfFilename = tempfile(fileext='.nc')
ncdfSolarFilename = tempfile(fileext='.nc')
# Build netCDF grids and over a defined time period.
# Only precip data is to be added to the netCDF files.
# This is because the URLs for the other variables are set to zero.
file.names = makeNetCDF_file(ncdfFilename=ncdfFilename,
ncdfSolarFilename=ncdfSolarFilename,
updateFrom=startDate, updateTo=endDate,
urlTmin=NA, urlTmax=NA, urlVprp=NA, urlSolarrad=NA)
# Load example catchment boundaries and remove all but the first.
# Note, this is done only to speed up the example runtime.
data("catchments")
catchments = catchments[1,]
# Extract daily precip. data (not Tmin, Tmax, VPD, ET).
# Note, the input "catchments" can also be a file to a ESRI shape file.
climateData = extractCatchmentData(ncdfFilename=file.names$ncdfFilename,
ncdfSolarFilename=file.names$ncdfSolarFilename,
extractFrom=startDate, extractTo=endDate,
getTmin = FALSE, getTmax = FALSE, getVprp = FALSE,
getSolarrad = FALSE, getET = FALSE,
catchments=catchments,
temporal.timestep = 'daily')
# Extract the daily catchment average data.
climateDataAvg = climateData$catchmentTemporal.mean
# Extract the daily catchment variance data.
climateDataVar = climateData$catchmentSpatial.var
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.