Collect USGS stream observations to evaluate streamflow simulation

Aubrey Dugger, James McCreight, Alyssa Hendricks

2017-05-01

Background

USGS streamflow observations are a primary source of hydrologic information and often used for validation and calibration of hydrlogic models. Recently, web services have been developed at NWIS and the dataRetrieval R package has emerged to make it easy to get USGS data into R.

This vignette demonstates some rwrfhydro tools to collect, store, and mainipulate USGS within a local database. This data is then used to evaluate model performance in estimating streamflow with and without overland, subsurface, and groundwater routing active. WRF-Hydro was run predict streamflow for Fourmile Creek at the Orodell USGS gage for the 2013 snowmelt period. NoahMP was used as the LSM and a 3-year spinup period was completed prior to the part of the simulation in the test case directory.

The fundamental layout of the local database is simply a directory containing:

Generally, there are two basic kinds of functions: “Get” and “Query”. Get functions use dataRetrieval functions to actually go out to NWIS and “get” data and metadata. Query functions query the local database. There are exceptions to this.

Setup

Load the rwrfhydro package.

library("rwrfhydro")

Create a directory and path to where you want your database to be built:

dbPath <- '~/wrfHydroTestCases/usgsDb/'

Discover gage locations, get data, and save to local database.

Perhaps you know the lat/lon of a gage (e.g. from frxst_pts_out.txt) and you need the HUC8. The within argument is taken to be in decimal degrees.

stnDf <- FindUsgsStns(stnLon=254.67374999999998408,
                      stnLat=40.018666670000001773,
                      within=.001)
str(stnDf)
## 'data.frame':    1 obs. of  13 variables:
##  $ stnLon       : num 255
##  $ stnLat       : num 40
##  $ within       : num 0.001
##  $ siteType     : chr "ST"
##  $ hasDataTypeCd: chr "iv"
##  $ agency_cd    : chr "USGS"
##  $ site_no      : chr "06727500"
##  $ station_nm   : chr "FOURMILE CREEK AT ORODELL, CO"
##  $ site_tp_cd   : chr "ST"
##  $ dec_lat_va   : num 40
##  $ dec_long_va  : num -105
##  $ colocated    : logi FALSE
##  $ queryTime    : POSIXct, format: "2017-05-01 16:43:05"

Now you have the site_no or station ID number, “06727500”. You could have also used FindUsgsStns() to reveal the gages in the HUC8, if you knew the HUC8 code. Because we organize the database by HUC8, we have a function to get HUC8 from station id. Then we get the above information for all locations in the HUC8.

huc8 <- GetSiteHuc(stnDf$site_no)
str(FindUsgsStns(huc=huc8))
## 'data.frame':    8 obs. of  11 variables:
##  $ huc8         : Factor w/ 1 level "10190005": 1 1 1 1 1 1 1 1
##  $ siteType     : chr  "ST" "ST" "ST" "ST" ...
##  $ hasDataTypeCd: chr  "iv" "iv" "iv" "iv" ...
##  $ agency_cd    : chr  "USGS" "USGS" "USGS" "USGS" ...
##  $ site_no      : chr  "06727410" "06730160" "06727500" "06730200" ...
##  $ station_nm   : chr  "FOURMILE CREEK AT LOGAN MILL ROAD NEAR CRISMAN, CO" "FOURMILE CANYON CREEK NEAR SUNSHINE, CO" "FOURMILE CREEK AT ORODELL, CO" "BOULDER CREEK AT NORTH 75TH ST. NEAR BOULDER, CO" ...
##  $ site_tp_cd   : chr  "ST" "ST" "ST" "ST" ...
##  $ dec_lat_va   : num  40 40.1 40 40.1 40.1 ...
##  $ dec_long_va  : num  -105 -105 -105 -105 -105 ...
##  $ colocated    : logi  FALSE FALSE FALSE FALSE FALSE FALSE ...
##  $ queryTime    : POSIXct, format: "2017-05-01 16:43:05" "2017-05-01 16:43:05" "2017-05-01 16:43:05" "2017-05-01 16:43:05" ...

FindUsgsStns is a wrapper on dataRetrieval::whatNWISsites which has been written to focus on instantaneous values. (It is worth noting the flexibility and generality of the underlying function.)

Now pull the data for this HUC8. Currently, this grabs all available products for the HUC. Note that the HUC data are organized by product code (e.g. 00060) then by data and meta in the returned list. (Also note that this command sometimes fails on the remote end and may need to be rerun.) In meta, siteInfo is the meta that we use in querying the local data base in commands shown below.

coData <- GetUsgsHucData(huc=huc8)  ## this can take a little while...
str(coData)
str(coData)
## List of 2
##  $ 00060:List of 2
##   ..$ data:'data.frame': 1252525 obs. of  6 variables:
##   .. ..$ agency_cd       : chr [1:1252525] "USGS" "USGS" "USGS" "USGS" ...
##   .. ..$ site_no         : chr [1:1252525] "06724970" "06724970" "06724970" "06724970" ...
##   .. ..$ dateTime        : POSIXct[1:1252525], format: "2014-03-06 07:00:00" "2014-03-06 07:10:00" "2014-03-06 07:20:00" "2014-03-06 07:30:00" ...
##   .. ..$ tz_cd           : chr [1:1252525] "UTC" "UTC" "UTC" "UTC" ...
##   .. ..$ X_00060_00011_cd: chr [1:1252525] "A" "A" "A" "A" ...
##   .. ..$ X_00060_00011   : num [1:1252525] 11 11 11 12 12 13 13 13 13 13 ...
##   ..$ meta:List of 3
##   .. ..$ siteInfo     :'data.frame': 8 obs. of  46 variables:
##   .. .. ..$ site_no              : chr [1:8] "06724970" "06725450" "06727410" "06727500" ...
##   .. .. ..$ agency_cd            : chr [1:8] "USGS" "USGS" "USGS" "USGS" ...
##   .. .. ..$ dec_lat_va           : num [1:8] 40.1 40.2 40 40 40.1 ...
##   .. .. ..$ dec_long_va          : num [1:8] -105 -105 -105 -105 -105 ...
##   .. .. ..$ huc_cd               : chr [1:8] "10190005" "10190005" "10190005" "10190005" ...
##   .. .. ..$ site_tp_cd           : chr [1:8] "ST" "ST" "ST" "ST" ...
##   .. .. ..$ state_cd             : chr [1:8] "08" "08" "08" "08" ...
##   .. .. ..$ station_nm           : chr [1:8] "LEFT HAND CREEK AT HOVER ROAD NEAR LONGMONT, CO" "ST. VRAIN CREEK BELOW LONGMONT, CO" "FOURMILE CREEK AT LOGAN MILL ROAD NEAR CRISMAN, CO" "FOURMILE CREEK AT ORODELL, CO" ...
##   .. .. ..$ tz_cd                : chr [1:8] "MST" "MST" "MST" "MST" ...
##   .. .. ..$ lat_va               : num [1:8] 400803 400927 400231 400107 400327 ...
##   .. .. ..$ long_va              : num [1:8] 1050751 1050055 1052154 1051934 1052056 ...
##   .. .. ..$ coord_meth_cd        : chr [1:8] "N" "N" "M" "N" ...
##   .. .. ..$ coord_acy_cd         : chr [1:8] "M" "5" "1" "5" ...
##   .. .. ..$ coord_datum_cd       : chr [1:8] "NAD83" "NAD83" "NAD83" "NAD83" ...
##   .. .. ..$ dec_coord_datum_cd   : chr [1:8] "NAD83" "NAD83" "NAD83" "NAD83" ...
##   .. .. ..$ district_cd          : chr [1:8] "08" "08" "08" "08" ...
##   .. .. ..$ county_cd            : chr [1:8] "013" "123" "013" "013" ...
##   .. .. ..$ country_cd           : chr [1:8] "US" "US" "US" "US" ...
##   .. .. ..$ land_net_ds          : chr [1:8] "  SWSWS16 T02N  R69W  6" "  NWNWS09 T002N R068W" "  NWNES20 T01N  R71W" "  NWSES27 T001N R071W" ...
##   .. .. ..$ map_nm               : chr [1:8] "HYGIENE, CO" "LONGMONT, CO" "BOULDER, CO" "BOULDER, CO" ...
##   .. .. ..$ map_scale_fc         : chr [1:8] "  24000" "  24000" "  24000" "  24000" ...
##   .. .. ..$ alt_va               : num [1:8] 4994 4852 6341 5750 6639 ...
##   .. .. ..$ alt_meth_cd          : chr [1:8] "N" "M" "N" "M" ...
##   .. .. ..$ alt_acy_va           : num [1:8] 4.3 10 4.3 15 4.3 10 10 4.3
##   .. .. ..$ alt_datum_cd         : chr [1:8] "NAVD88" "NGVD29" "NAVD88" "NGVD29" ...
##   .. .. ..$ basin_cd             : chr [1:8] "06" "06" "" "06" ...
##   .. .. ..$ topo_cd              : chr [1:8] "" "" "C" "" ...
##   .. .. ..$ instruments_cd       : chr [1:8] "YNNNYNNNNNNNNNNNNNNNNNNNNNNNNN" "NNNNYNNNNNNNNNNNNNNNNNNNNNNNNN" "NNNNYNYNNYNNNNNNNNNNNNNNNNNNNN" "NNNNYNYNNYNNNNNNNNNNNNNNNNNNNN" ...
##   .. .. ..$ construction_dt      : chr [1:8] "" "" "" "" ...
##   .. .. ..$ inventory_dt         : chr [1:8] "20140304" "" "20110317" "20110323" ...
##   .. .. ..$ drain_area_va        : num [1:8] 71.6 423 19.2 24.2 1.8 307 447 NA
##   .. .. ..$ contrib_drain_area_va: num [1:8] NA NA NA NA NA NA NA NA
##   .. .. ..$ local_time_fg        : chr [1:8] "Y" "Y" "Y" "Y" ...
##   .. .. ..$ reliability_cd       : chr [1:8] "U" "" "" "" ...
##   .. .. ..$ gw_file_cd           : chr [1:8] " Y" "NYNNNNNN" " Y" "NYNNNNNN" ...
##   .. .. ..$ nat_aqfr_cd          : chr [1:8] "" "" "" "" ...
##   .. .. ..$ aqfr_cd              : chr [1:8] "" "" "" "" ...
##   .. .. ..$ aqfr_type_cd         : chr [1:8] "" "" "" "" ...
##   .. .. ..$ well_depth_va        : num [1:8] NA NA NA NA NA NA NA NA
##   .. .. ..$ hole_depth_va        : num [1:8] NA NA NA NA NA NA NA NA
##   .. .. ..$ depth_src_cd         : chr [1:8] "" "" "" "" ...
##   .. .. ..$ project_no           : chr [1:8] "00001200" "" "000012000" "000012000" ...
##   .. .. ..$ timeZoneOffset       : chr [1:8] "-07:00" "-07:00" "-07:00" "-07:00" ...
##   .. .. ..$ srs                  : chr [1:8] "EPSG:4326" "EPSG:4326" "EPSG:4326" "EPSG:4326" ...
##   .. .. ..$ startTime            : POSIXct[1:8], format: "2014-03-06 07:00:00" "2007-10-01 06:00:00" "2011-04-01 07:00:00" "2011-04-01 06:00:00" ...
##   .. .. ..$ endTime              : POSIXct[1:8], format: "2015-04-17 04:10:00" "2013-09-12 05:45:00" "2013-09-12 05:55:00" "2015-04-17 04:10:00" ...
##   .. ..$ variableInfo :'data.frame': 1 obs. of  6 variables:
##   .. .. ..$ parameterCd   : chr "00060"
##   .. .. ..$ parameter_nm  : chr "Streamflow, ft&#179;/s"
##   .. .. ..$ parameter_desc: chr "Discharge, cubic feet per second"
##   .. .. ..$ valueType     : chr "Derived Value"
##   .. .. ..$ param_units   : chr "ft3/s"
##   .. .. ..$ noDataValue   : logi NA
##   .. ..$ statisticInfo:'data.frame': 1 obs. of  2 variables:
##   .. .. ..$ statisticName: chr ""
##   .. .. ..$ statisticCd  : chr "00011"
##  $ 00065:List of 2
##   ..$ data:'data.frame': 64142 obs. of  6 variables:
##   .. ..$ agency_cd       : chr [1:64142] "USGS" "USGS" "USGS" "USGS" ...
##   .. ..$ site_no         : chr [1:64142] "06724970" "06724970" "06724970" "06724970" ...
##   .. ..$ dateTime        : POSIXct[1:64142], format: "2014-12-17 07:00:00" "2014-12-17 07:10:00" "2014-12-17 07:20:00" "2014-12-17 07:30:00" ...
##   .. ..$ tz_cd           : chr [1:64142] "UTC" "UTC" "UTC" "UTC" ...
##   .. ..$ X_00065_00011_cd: chr [1:64142] "P" "P" "P" "P" ...
##   .. ..$ X_00065_00011   : num [1:64142] 6.48 6.49 6.48 6.49 6.49 6.49 6.48 6.48 6.48 6.48 ...
##   ..$ meta:List of 3
##   .. ..$ siteInfo     :'data.frame': 6 obs. of  46 variables:
##   .. .. ..$ site_no              : chr [1:6] "06724970" "06727500" "06730160" "06730200" ...
##   .. .. ..$ agency_cd            : chr [1:6] "USGS" "USGS" "USGS" "USGS" ...
##   .. .. ..$ dec_lat_va           : num [1:6] 40.1 40 40.1 40.1 40.1 ...
##   .. .. ..$ dec_long_va          : num [1:6] -105 -105 -105 -105 -105 ...
##   .. .. ..$ huc_cd               : chr [1:6] "10190005" "10190005" "10190005" "10190005" ...
##   .. .. ..$ site_tp_cd           : chr [1:6] "ST" "ST" "ST" "ST" ...
##   .. .. ..$ state_cd             : chr [1:6] "08" "08" "08" "08" ...
##   .. .. ..$ station_nm           : chr [1:6] "LEFT HAND CREEK AT HOVER ROAD NEAR LONGMONT, CO" "FOURMILE CREEK AT ORODELL, CO" "FOURMILE CANYON CREEK NEAR SUNSHINE, CO" "BOULDER CREEK AT NORTH 75TH ST. NEAR BOULDER, CO" ...
##   .. .. ..$ tz_cd                : chr [1:6] "MST" "MST" "MST" "MST" ...
##   .. .. ..$ lat_va               : num [1:6] 400803 400107 400327 400306 400820 ...
##   .. .. ..$ long_va              : num [1:6] 1050751 1051934 1052056 1051042 1050113 ...
##   .. .. ..$ coord_meth_cd        : chr [1:6] "N" "N" "G" "M" ...
##   .. .. ..$ coord_acy_cd         : chr [1:6] "M" "5" "S" "M" ...
##   .. .. ..$ coord_datum_cd       : chr [1:6] "NAD83" "NAD83" "NAD83" "NAD27" ...
##   .. .. ..$ dec_coord_datum_cd   : chr [1:6] "NAD83" "NAD83" "NAD83" "NAD83" ...
##   .. .. ..$ district_cd          : chr [1:6] "08" "08" "08" "08" ...
##   .. .. ..$ county_cd            : chr [1:6] "013" "013" "013" "013" ...
##   .. .. ..$ country_cd           : chr [1:6] "US" "US" "US" "US" ...
##   .. .. ..$ land_net_ds          : chr [1:6] "  SWSWS16 T02N  R69W  6" "  NWSES27 T001N R071W" "      S09 T01N  R71W  6" "  SENWS13 T001N R070W" ...
##   .. .. ..$ map_nm               : chr [1:6] "HYGIENE, CO" "BOULDER, CO" "BOULDER, CO" "" ...
##   .. .. ..$ map_scale_fc         : chr [1:6] "  24000" "  24000" "  24000" "" ...
##   .. .. ..$ alt_va               : num [1:6] 4994 5750 6639 5106 4860 ...
##   .. .. ..$ alt_meth_cd          : chr [1:6] "N" "M" "N" "M" ...
##   .. .. ..$ alt_acy_va           : num [1:6] 4.3 15 4.3 10 10 4.3
##   .. .. ..$ alt_datum_cd         : chr [1:6] "NAVD88" "NGVD29" "NAVD88" "NGVD29" ...
##   .. .. ..$ basin_cd             : chr [1:6] "06" "06" "06" "" ...
##   .. .. ..$ topo_cd              : chr [1:6] "" "" "" "" ...
##   .. .. ..$ instruments_cd       : chr [1:6] "YNNNYNNNNNNNNNNNNNNNNNNNNNNNNN" "NNNNYNYNNYNNNNNNNNNNNNNNNNNNNN" "NNNNYNNNNYNNNNNNNNNNNNNNNNNNNN" "YNNNYNNNNNNNNNNNNNNNNNNNNNNNNN" ...
##   .. .. ..$ construction_dt      : chr [1:6] "" "" "20120416" "" ...
##   .. .. ..$ inventory_dt         : chr [1:6] "20140304" "20110323" "" "" ...
##   .. .. ..$ drain_area_va        : num [1:6] 71.6 24.2 1.8 307 447 NA
##   .. .. ..$ contrib_drain_area_va: num [1:6] NA NA NA NA NA NA
##   .. .. ..$ local_time_fg        : chr [1:6] "Y" "Y" "Y" "Y" ...
##   .. .. ..$ reliability_cd       : chr [1:6] "U" "" "" "" ...
##   .. .. ..$ gw_file_cd           : chr [1:6] " Y" "NYNNNNNN" " Y" "NYNNNNNN" ...
##   .. .. ..$ nat_aqfr_cd          : chr [1:6] "" "" "" "" ...
##   .. .. ..$ aqfr_cd              : chr [1:6] "" "" "" "" ...
##   .. .. ..$ aqfr_type_cd         : chr [1:6] "" "" "" "" ...
##   .. .. ..$ well_depth_va        : num [1:6] NA NA NA NA NA NA
##   .. .. ..$ hole_depth_va        : num [1:6] NA NA NA NA NA NA
##   .. .. ..$ depth_src_cd         : chr [1:6] "" "" "" "" ...
##   .. .. ..$ project_no           : chr [1:6] "00001200" "000012000" "000012000" "" ...
##   .. .. ..$ timeZoneOffset       : chr [1:6] "-07:00" "-07:00" "-07:00" "-07:00" ...
##   .. .. ..$ srs                  : chr [1:6] "EPSG:4326" "EPSG:4326" "EPSG:4326" "EPSG:4326" ...
##   .. .. ..$ startTime            : POSIXct[1:6], format: "2014-12-17 07:00:00" "2015-03-27 16:15:00" "2015-03-27 19:35:00" "2014-12-17 07:00:00" ...
##   .. .. ..$ endTime              : POSIXct[1:6], format: "2015-04-17 04:10:00" "2015-04-17 04:10:00" "2015-04-17 04:20:00" "2015-04-17 04:15:00" ...
##   .. ..$ variableInfo :'data.frame': 1 obs. of  6 variables:
##   .. .. ..$ parameterCd   : chr "00065"
##   .. .. ..$ parameter_nm  : chr "Gage height, ft"
##   .. .. ..$ parameter_desc: chr "Gage height, feet"
##   .. .. ..$ valueType     : chr "Derived Value"
##   .. .. ..$ param_units   : chr "ft"
##   .. .. ..$ noDataValue   : logi NA
##   .. ..$ statisticInfo:'data.frame': 1 obs. of  2 variables:
##   .. .. ..$ statisticName: chr ""
##   .. .. ..$ statisticCd  : chr "00011"
##  - attr(*, "split_type")= chr "data.frame"
##  - attr(*, "split_labels")='data.frame': 2 obs. of  1 variable:
##   ..$ parm_cd: chr [1:2] "00060" "00065"

Now save this data to the local database. Note that this also could have been achieved by specifying the outPath argument to GetUsgsHucData.

coFiles <- SaveHucData(coData, outPath=dbPath)
coFiles

Query the local data

Now we work entirely locally, having grabbed the data of interest. For all HUC8 and products, any of the siteInfo metadata can be retrieved from the local DB. Note that the same site is repeated for multiple products.

QuerySiteInfo(c('station_nm','site_no','dec_lat_va','dec_long_va'), path=dbPath)
station_nm site_no dec_lat_va dec_long_va product HUC8
LEFT HAND CREEK AT HOVER ROAD NEAR LONGMONT, CO 06724970 40.13428 -105.13082 00060 10190005
ST. VRAIN CREEK BELOW LONGMONT, CO 06725450 40.15742 -105.01539 00060 10190005
FOURMILE CREEK AT LOGAN MILL ROAD NEAR CRISMAN, CO 06727410 40.04203 -105.36492 00060 10190005
FOURMILE CREEK AT ORODELL, CO 06727500 40.01867 -105.32625 00060 10190005
FOURMILE CANYON CREEK NEAR SUNSHINE, CO 06730160 40.05761 -105.34878 00060 10190005
BOULDER CREEK AT NORTH 75TH ST. NEAR BOULDER, CO 06730200 40.05165 -105.17888 00060 10190005
BOULDER CREEK AT MOUTH NEAR LONGMONT, CO 06730500 40.13878 -105.02022 00060 10190005
ST VRAIN CR BLW BOULDER CR AT HWY 119 NR LONGMONT 06730525 40.16035 -105.00783 00060 10190005
LEFT HAND CREEK AT HOVER ROAD NEAR LONGMONT, CO 06724970 40.13428 -105.13082 00065 10190005
FOURMILE CREEK AT ORODELL, CO 06727500 40.01867 -105.32625 00065 10190005
FOURMILE CANYON CREEK NEAR SUNSHINE, CO 06730160 40.05761 -105.34878 00065 10190005
BOULDER CREEK AT NORTH 75TH ST. NEAR BOULDER, CO 06730200 40.05165 -105.17888 00065 10190005
BOULDER CREEK AT MOUTH NEAR LONGMONT, CO 06730500 40.13878 -105.02022 00065 10190005
ST VRAIN CR BLW BOULDER CR AT HWY 119 NR LONGMONT 06730525 40.16035 -105.00783 00065 10190005
TOMBIGBEE R AT COFFEEVILLE L&D NR COFFEEVILLE, AL. 02469761 31.75849 -88.12918 00060 03160203
SATILPA CREEK NEAR COFFEEVILLE AL 02469800 31.74432 -88.02251 00060 03160203
BASSETT CREEK AT US HIGHWAY 43 NR THOMASVILLE, AL. 02470072 31.86404 -87.74722 00060 03160203
TOMBIGBEE R AT COFFEEVILLE L&D NR COFFEEVILLE, AL. 02469761 31.75849 -88.12918 00065 03160203
TOMBIGBEE R BL COFFEEVILLE L&D NEAR COFFEEVILLE 02469762 31.75710 -88.12501 00065 03160203
SATILPA CREEK NEAR COFFEEVILLE AL 02469800 31.74432 -88.02251 00065 03160203
TOMBIGBEE RIVER AT STEAMPLANT NR LEROY, AL 02470050 31.48600 -87.90889 00065 03160203
BASSETT CREEK AT US HIGHWAY 43 NR THOMASVILLE, AL. 02470072 31.86404 -87.74722 00065 03160203
TOMBIGBEE R AT COFFEEVILLE L&D NR COFFEEVILLE, AL. 02469761 31.75849 -88.12918 45592 03160203

Say you just want Orodell and you want your code to be readable: translate the name to the code with QuerySiteName (which translates both ways).

dataOrodell <- QuerySiteData(QuerySiteName("FOURMILE CREEK AT ORODELL, CO", path=dbPath),
                             product='00060', path=dbPath)
str(dataOrodell)
## 'data.frame':    208824 obs. of  5 variables:
##  $ site_no : chr  "06727500" "06727500" "06727500" "06727500" ...
##  $ dateTime: POSIXct, format: "2011-04-01 06:00:00" "2011-04-01 06:05:00" "2011-04-01 06:10:00" "2011-04-01 06:15:00" ...
##  $ code    : chr  "A" "A" "A" "A" ...
##  $ value   : num  1.2 1.2 1.2 1.2 1.2 1.2 1.2 1.2 1.2 1.2 ...
##  $ variable: chr  "X_00060_00011" "X_00060_00011" "X_00060_00011" "X_00060_00011" ...
##  - attr(*, "siteInfo")='data.frame': 1 obs. of  46 variables:
##   ..$ site_no              : chr "06727500"
##   ..$ agency_cd            : chr "USGS"
##   ..$ dec_lat_va           : num 40
##   ..$ dec_long_va          : num -105
##   ..$ huc_cd               : chr "10190005"
##   ..$ site_tp_cd           : chr "ST"
##   ..$ state_cd             : chr "08"
##   ..$ station_nm           : chr "FOURMILE CREEK AT ORODELL, CO"
##   ..$ tz_cd                : chr "MST"
##   ..$ lat_va               : num 4e+05
##   ..$ long_va              : num 1051934
##   ..$ coord_meth_cd        : chr "N"
##   ..$ coord_acy_cd         : chr "5"
##   ..$ coord_datum_cd       : chr "NAD83"
##   ..$ dec_coord_datum_cd   : chr "NAD83"
##   ..$ district_cd          : chr "08"
##   ..$ county_cd            : chr "013"
##   ..$ country_cd           : chr "US"
##   ..$ land_net_ds          : chr "  NWSES27 T001N R071W"
##   ..$ map_nm               : chr "BOULDER, CO"
##   ..$ map_scale_fc         : chr "  24000"
##   ..$ alt_va               : num 5750
##   ..$ alt_meth_cd          : chr "M"
##   ..$ alt_acy_va           : num 15
##   ..$ alt_datum_cd         : chr "NGVD29"
##   ..$ basin_cd             : chr "06"
##   ..$ topo_cd              : chr ""
##   ..$ instruments_cd       : chr "NNNNYNYNNYNNNNNNNNNNNNNNNNNNNN"
##   ..$ construction_dt      : chr ""
##   ..$ inventory_dt         : chr "20110323"
##   ..$ drain_area_va        : num 24.2
##   ..$ contrib_drain_area_va: num NA
##   ..$ local_time_fg        : chr "Y"
##   ..$ reliability_cd       : chr ""
##   ..$ gw_file_cd           : chr "NYNNNNNN"
##   ..$ nat_aqfr_cd          : chr ""
##   ..$ aqfr_cd              : chr ""
##   ..$ aqfr_type_cd         : chr ""
##   ..$ well_depth_va        : num NA
##   ..$ hole_depth_va        : num NA
##   ..$ depth_src_cd         : chr ""
##   ..$ project_no           : chr "000012000"
##   ..$ timeZoneOffset       : chr "-07:00"
##   ..$ srs                  : chr "EPSG:4326"
##   ..$ startTime            : POSIXct, format: "2011-04-01 06:00:00"
##   ..$ endTime              : POSIXct, format: "2015-04-17 04:10:00"
##  - attr(*, "variableInfo")='data.frame': 1 obs. of  6 variables:
##   ..$ parameterCd   : chr "00060"
##   ..$ parameter_nm  : chr "Streamflow, ft&#179;/s"
##   ..$ parameter_desc: chr "Discharge, cubic feet per second"
##   ..$ valueType     : chr "Derived Value"
##   ..$ param_units   : chr "ft3/s"
##   ..$ noDataValue   : logi NA
##  - attr(*, "statisticInfo")='data.frame':    1 obs. of  2 variables:
##   ..$ statisticName: chr ""
##   ..$ statisticCd  : chr "00011"

Now make it “pretty”. The main difference here is meaningful column names and identification of variables and codes in the attributes. We have defined “prettyUsgs” as an S3 class.

prettyOrodell <- PrettyUsgs(dataOrodell, metric=TRUE)
str(prettyOrodell)
## Classes 'prettyUsgs' and 'data.frame':   207496 obs. of  5 variables:
##  $ site_no : chr  "06727500" "06727500" "06727500" "06727500" ...
##  $ POSIXct : POSIXct, format: "2011-04-01 06:00:00" "2011-04-01 06:05:00" "2011-04-01 06:10:00" "2011-04-01 06:15:00" ...
##  $ code    : chr  "A" "A" "A" "A" ...
##  $ value   : num  0.034 0.034 0.034 0.034 0.034 ...
##  $ variable: chr  "Discharge (cms)" "Discharge (cms)" "Discharge (cms)" "Discharge (cms)" ...
##  - attr(*, "siteInfo")='data.frame': 1 obs. of  46 variables:
##   ..$ site_no              : chr "06727500"
##   ..$ agency_cd            : chr "USGS"
##   ..$ dec_lat_va           : num 40
##   ..$ dec_long_va          : num -105
##   ..$ huc_cd               : chr "10190005"
##   ..$ site_tp_cd           : chr "ST"
##   ..$ state_cd             : chr "08"
##   ..$ station_nm           : chr "FOURMILE CREEK AT ORODELL, CO"
##   ..$ tz_cd                : chr "MST"
##   ..$ lat_va               : num 4e+05
##   ..$ long_va              : num 1051934
##   ..$ coord_meth_cd        : chr "N"
##   ..$ coord_acy_cd         : chr "5"
##   ..$ coord_datum_cd       : chr "NAD83"
##   ..$ dec_coord_datum_cd   : chr "NAD83"
##   ..$ district_cd          : chr "08"
##   ..$ county_cd            : chr "013"
##   ..$ country_cd           : chr "US"
##   ..$ land_net_ds          : chr "  NWSES27 T001N R071W"
##   ..$ map_nm               : chr "BOULDER, CO"
##   ..$ map_scale_fc         : chr "  24000"
##   ..$ alt_va               : num 5750
##   ..$ alt_meth_cd          : chr "M"
##   ..$ alt_acy_va           : num 15
##   ..$ alt_datum_cd         : chr "NGVD29"
##   ..$ basin_cd             : chr "06"
##   ..$ topo_cd              : chr ""
##   ..$ instruments_cd       : chr "NNNNYNYNNYNNNNNNNNNNNNNNNNNNNN"
##   ..$ construction_dt      : chr ""
##   ..$ inventory_dt         : chr "20110323"
##   ..$ drain_area_va        : num 24.2
##   ..$ contrib_drain_area_va: num NA
##   ..$ local_time_fg        : chr "Y"
##   ..$ reliability_cd       : chr ""
##   ..$ gw_file_cd           : chr "NYNNNNNN"
##   ..$ nat_aqfr_cd          : chr ""
##   ..$ aqfr_cd              : chr ""
##   ..$ aqfr_type_cd         : chr ""
##   ..$ well_depth_va        : num NA
##   ..$ hole_depth_va        : num NA
##   ..$ depth_src_cd         : chr ""
##   ..$ project_no           : chr "000012000"
##   ..$ timeZoneOffset       : chr "-07:00"
##   ..$ srs                  : chr "EPSG:4326"
##   ..$ startTime            : POSIXct, format: "2011-04-01 06:00:00"
##   ..$ endTime              : POSIXct, format: "2015-04-17 04:10:00"
class(prettyOrodell)
## [1] "prettyUsgs" "data.frame"

Plot the “pretty”" data.

oroPlot <- PlotPrettyUsgs(prettyOrodell)

Import modelled and observed datasets

Set a data path to the Fourmile Creek test case.

dataPath <- '~/wrfHydroTestCases/Fourmile_Creek_testcase_v2.0/'

Model 1: Only channel routing turned on (hourly model run).

modStrd.chrt.fc <- ReadFrxstPts(paste0(dataPath, '/run.ChannelRouting/frxst_pts_out.txt'))

Model 2: All WRF-Hydro routing options turned on (hourly model run).

modStrd.allrt.fc <- ReadFrxstPts(paste0(dataPath, '/run.FullRouting/frxst_pts_out.txt'))

USGS gage observed data at 5-minute intervals. Find the nearest gage to the forecast point in the above files. (The following approach can be used when multiple forecast points are output by the model. The ddply function essentially performs a loop over unique st_id passing the corresponding subset of the data frame modStrd.chrt.fc to the function which just returns the first row’s lon and lat.)

library(plyr)
fcLocation <- ddply(modStrd.chrt.fc, .(st_id), function(df) df[1,c('st_lon','st_lat')])
near <- FindUsgsStns(stnLon=c(fcLocation$st_lon), stnLat=c(fcLocation$st_lat), within=.005)
near
stnLon stnLat within siteType hasDataTypeCd agency_cd site_no station_nm site_tp_cd dec_lat_va dec_long_va colocated queryTime
-105.3278 40.01932 0.005 ST iv USGS 06727500 FOURMILE CREEK AT ORODELL, CO ST 40.01867 -105.3263 FALSE 2017-05-01 16:43:19

The following checks to see if the data are local and gets them if not. The ‘00060’ product (streamflow) is returned and then made pretty, including conversion to metric.

dbPath <- '~/wrfHydroTestCases/usgsDb/'
obsStr5min.fc <- 
  PrettyUsgs(QueryHaveSite(near$site_no, path=dbPath, ret='00060', get=TRUE))
obsStr5min.fc <- plyr::rename(obsStr5min.fc, c(value='q_cms'))

Plot hydrographs

Compare hydrographs for the full model run.

PlotFluxCompare(obsStr5min.fc, "q_cms", modStrd.chrt.fc, "q_cms", strDf.mod2=modStrd.allrt.fc, 
     strCol.mod2="q_cms", labelObs="Observed Fourmile Creek at Orodell", 
     labelMod1="Channel Routing Only", labelMod2="All Routing", 
     title="Streamflow: Fourmile Creek")

Now limit the plot to the peak May flow period only. The reported stats are updated to the new time period. (Note that the R warning is innocuous because the subset adjusts for timezone, so it’s ok that the timezones don’t match.)

PlotFluxCompare(obsStr5min.fc, "q_cms", modStrd.chrt.fc, "q_cms", strDf.mod2=modStrd.allrt.fc, 
     strCol.mod2="q_cms", labelObs="Observed Fourmile Creek at Orodell", 
     labelMod1="Channel Routing Only", labelMod2="All Routing", 
     title="Streamflow: Fourmile Creek", 
     stdate=as.POSIXct("2013-05-01 00:00:00", format="%Y-%m-%d %H:%M:%S", tz="UTC"), 
     enddate=as.POSIXct("2013-05-31 00:00:00", format="%Y-%m-%d %H:%M:%S", tz="UTC"))

Review flow duration curves

NOTE: You generally evaluate flow duration curves and staistics over much longer time periods (e.g., multiple years) than what we demo here. To make the test case more portable, we are only evaluating once-a-day model output over 5 months.

Calculate percent exceedances for flow duration curves. Note that we need to subset the observations to match our model run output times, and vice versa.

obsStr5min.comp.fc <- CalcFdc(subset(obsStr5min.fc, POSIXct %in% c(modStrd.chrt.fc$POSIXct)))
modStrd.chrt.comp.fc <- CalcFdc(subset(modStrd.chrt.fc, POSIXct %in% c(obsStr5min.comp.fc$POSIXct)))
modStrd.allrt.comp.fc <- CalcFdc(subset(modStrd.allrt.fc, POSIXct %in% c(obsStr5min.comp.fc$POSIXct)))

Compare how the models are doing predicting flow values that will be exceeded 20% of the time. First, calculate the fitted spline functions.

fdc.obsStr5min.comp.fc <- CalcFdcSpline(obsStr5min.comp.fc)
fdc.modStrd.chrt.comp.fc <- CalcFdcSpline(modStrd.chrt.comp.fc)
fdc.modStrd.allrt.comp.fc <- CalcFdcSpline(modStrd.allrt.comp.fc)

Then, evaluate at the 20% exceedance percentage (high flows).

fdc.obsStr5min.comp.fc(0.2)
fdc.modStrd.chrt.comp.fc(0.2)
fdc.modStrd.allrt.comp.fc(0.2)
## [1] 0.7706495
## [1] 0.8170807
## [1] 0.7931814

Now try the 80% exceedance percentage (low flows).

fdc.obsStr5min.comp.fc(0.8)
fdc.modStrd.chrt.comp.fc(0.8)
fdc.modStrd.allrt.comp.fc(0.8)
## [1] 0.08530777
## [1] 0.1557625
## [1] 0.09120762

Plot flow duration curves for a more complete picture. This tool will do the date matching for us, so no need to subset the datasets.

PlotFdcCompare(obsStr5min.fc, "q_cms", modStrd.chrt.fc, "q_cms", strDf.mod2=modStrd.allrt.fc, 
     strCol.mod2="q_cms", labelObs="Observed Fourmile Creek", 
     labelMod1="Channel Routing Only", labelMod2="All Routing")
## [1] "Combined min flow for y-axis (capped at 0.001): 0.03114853126"
## [1] "Combined max flow for y-axis: 2.5201993474"

Review model performance statistics

Calculate model performance stats (special formatting comands hidden). Again, this tool does the date matching for us.

CalcModPerf(modStrd.chrt.fc, obsStr5min.fc)
  n nse nselog cor rmse rmsenorm bias msd mae errcom errmaxt errfdc
units count unitless unitless unitless flux units unitless percent flux units flux units hours|days hours|days flux units
t 232 0.82 0.73 0.93 0.2067 8.31 19.2 0.0846 0.1305 NA NA 0.08
daily 116 0.82 0.73 0.93 0.2067 8.31 19.2 0.0846 0.1305 0 0 0.08
monthly 5 0.92 0.7 0.99 0.1044 10.07 21.5 0.0768 0.0806 -0.2 -0.2 NA
yearly 1 NA NA NA NA NA 19.2 0.0846 0.0846 3 -1 NA
max10 24 -0.04 -0.08 0.5 0.4806 32.64 15 0.235 0.417 NA NA NA
min10 24 -82.69 -24.6 0.53 0.078 250.56 152.6 0.076 0.076 NA NA NA
CalcModPerf(modStrd.allrt.fc, obsStr5min.fc)
  n nse nselog cor rmse rmsenorm bias msd mae errcom errmaxt errfdc
units count unitless unitless unitless flux units unitless percent flux units flux units hours|days hours|days flux units
t 116 0.75 0.87 0.9 0.2462 9.89 9.1 0.0399 0.1402 NA NA 0.03
daily 116 0.75 0.87 0.9 0.2462 9.89 9.1 0.0399 0.1402 0 0 0.03
monthly 5 0.95 0.95 0.99 0.0867 8.36 7 0.0251 0.0645 -0.6 0 NA
yearly 1 NA NA NA NA NA 9.1 0.0399 0.0399 3 4 NA
max10 12 -0.63 -0.51 0.11 0.6034 40.98 14.4 0.2264 0.5365 NA NA NA
min10 12 -13.12 -6.56 -0.08 0.0321 102.9 49.8 0.0248 0.0252 NA NA NA

Help on CalcModPerf gives details on the individual statistics returned.

help(CalcModPerf)
## Using development documentation for CalcModPerf
## Using development documentation for CalcModPerf



Calculate flow duration curve performance statistics.

CalcFdcPerf(modStrd.allrt.fc, obsStr5min.fc)
p.exceed q.mod q.obs q.err q.perr
0.1 1.43 1.04 0.39 37.5
0.2 0.79 0.77 0.02 2.6
0.3 0.49 0.56 -0.07 -12.5
0.4 0.33 0.36 -0.03 -8.3
0.5 0.25 0.24 0.01 4.2
0.6 0.17 0.14 0.03 21.4
0.7 0.13 0.12 0.01 8.3
0.8 0.09 0.09 0 0
0.9 0.07 0.06 0.01 16.7

Again, help on CalcFdcPerf gives details on the individual statistics returned.

help(CalcFdcPerf)
## Using development documentation for CalcFdcPerf
## Using development documentation for CalcFdcPerf

GagesII Attributes

We’ve imported the gages-II atributes to be directly available in R.

head(gages2Attr)
STAID STANAME CLASS AGGECOREGI DRAIN_SQKM HUC02 LAT_GAGE LNG_GAGE STATE HCDN_2009 ACTIVE09 FLYRS1900 FLYRS1950 FLYRS1990
01011000 Allagash River near Allagash, Maine Non-ref NorthEast 3186.8440 01 47.06961 -69.07954 ME yes 78 60 20
01013500 Fish River near Fort Kent, Maine Ref NorthEast 2252.6960 01 47.23739 -68.58264 ME yes yes 85 60 20
01015800 Aroostook River near Masardis, Maine Non-ref NorthEast 2313.7550 01 46.52300 -68.37176 ME yes 51 51 19
01016500 MACHIAS RIVER NEAR ASHLAND, ME Non-ref NorthEast 847.7721 01 46.62831 -68.43479 ME 32 32 0
01017000 Aroostook River at Washburn, Maine Non-ref NorthEast 4278.9070 01 46.77729 -68.15719 ME yes 79 60 20
01017060 Hardwood Brook below Glidden Brk nr Caribou, Maine Non-ref NorthEast 14.8644 01 46.78361 -67.98917 ME yes 0 0 0
?gages2Attr
## Using development documentation for gages2Attr