Skip to contents

Creates a data.table from a rotated pole netcdf (as usually found in RCMs), which includes values and date. Useful for extracting e.g. the series for a station. Requires that dimension variables in netcdf file contain rlon and rlat, and that it contains daily data.

Usage

rotpole_nc_point_to_dt(
  filename,
  variable,
  point_lon,
  point_lat,
  interpolate_to_standard_calendar = FALSE,
  verbose = FALSE,
  add_grid_coord = FALSE
)

Arguments

filename

Complete path to .nc file.

variable

Name of the variable to extract from filename (character).

point_lon

Numeric longitude of the point to extract (decimal degrees).

point_lat

Numeric latitude of the point to extract (decimal degrees).

interpolate_to_standard_calendar

Boolean, if TRUE will use map_non_standard_calendar to interpolate values to a standard calendar.

verbose

Boolean, if TRUE, will print more information.

add_grid_coord

Boolean, if TRUE, will add columns to the result which give the longitude and latitude of the underlying grid.

Value

A data.table with two columns: the dates in date, and the values in a variable named after input variable. The date column is of class Date, unless the .nc file has a non-standard calendar (360, noleap) and interpolate_to_standard_calendar is set to FALSE, in which it will be character. If add_grid_coord is set to TRUE, then two more columns named grid_lon and grid_lat.

Details

Calculates the euclidean distance, and takes the grid cell with minimal distance to point_lon and point_lat. Requires that the .nc file contains variables lon[rlon, rlat] and lat[rlon, rlat].

Examples

# example data from EURO-CORDEX (cropped for size)

# standard calendar
fn1 <- system.file("extdata", "test1.nc", package = "eurocordexr")
dt1 <- rotpole_nc_point_to_dt(
  filename = fn1,
  variable = "tasmin",
  point_lon = 11.31,
  point_lat = 46.5,
  verbose = TRUE
)
#> Succesfully opened file: /home/runner/work/_temp/Library/eurocordexr/extdata/test1.nc 
#> Point longitude =  11.31  ## Closest grid cell =  11.30084 
#> Point latitude =  46.5  ## Closest grid cell =  46.52603 
#> Euclidean distance in degrees =  0.02759755 

# non-standard calendar (360)
fn2 <- system.file("extdata", "test2.nc", package = "eurocordexr")

# read as is
dt2 <- rotpole_nc_point_to_dt(fn2, "tasmin", 11.31, 46.5)
str(dt2) # chr date
#> Classes ‘data.table’ and 'data.frame':	390 obs. of  2 variables:
#>  $ date  : chr  "1949-12-01" "1949-12-02" "1949-12-03" "1949-12-04" ...
#>  $ tasmin: num  264 263 262 270 275 ...
#>  - attr(*, ".internal.selfref")=<externalptr> 
dt2[86:94, ] # e.g. 30th of February in 360 calendar
#>          date   tasmin
#>        <char>    <num>
#> 1: 1950-02-26 273.2633
#> 2: 1950-02-27 272.6217
#> 3: 1950-02-28 266.5659
#> 4: 1950-02-29 262.9583
#> 5: 1950-02-30 264.9355
#> 6: 1950-03-01 268.0822
#> 7: 1950-03-02 266.9084
#> 8: 1950-03-03 267.9066
#> 9: 1950-03-04 266.9478

# interpolate to standard
dt3 <- rotpole_nc_point_to_dt(fn2, "tasmin", 11.31, 46.5,
                              interpolate_to_standard_calendar = TRUE)
str(dt3) # class Date
#> Classes ‘data.table’ and 'data.frame':	395 obs. of  2 variables:
#>  $ date  : Date, format: "1949-12-01" "1949-12-02" ...
#>  $ tasmin: num  264 263 262 270 275 ...
#>  - attr(*, ".internal.selfref")=<externalptr> 
dt3[86:94, ] # standard calender
#>          date   tasmin
#>        <Date>    <num>
#> 1: 1950-02-24 271.9780
#> 2: 1950-02-25 273.2633
#> 3: 1950-02-26 272.6217
#> 4: 1950-02-27 266.5659
#> 5: 1950-02-28 262.9583
#> 6: 1950-03-01 264.9355
#> 7: 1950-03-02 268.0822
#> 8: 1950-03-03 266.9084
#> 9: 1950-03-04 267.9066