cmsaf.cat {cmsafops} | R Documentation |
Concatenate datasets of several NetCDF input files.
Description
This function concatenates datasets of an arbitrary number of input files. All input files have to have the same structure with the same variable and different timesteps.
Usage
cmsaf.cat(var, infiles, outfile, nc34 = 4, overwrite = FALSE, verbose = FALSE)
Arguments
var |
Name of NetCDF variable (character). |
infiles |
Vector with filenames of input NetCDF files. The file names may include the directory (character). |
outfile |
Filename of output NetCDF file. This may include the directory (character). |
nc34 |
NetCDF version of output file. If |
overwrite |
logical; should existing output file be overwritten? |
verbose |
logical; if TRUE, progress messages are shown |
Value
A NetCDF file including the merged time series is written. The resulting file uses the meta data of the first input file.
Examples
## Create an example NetCDF file with a similar structure as used by CM
## SAF. The file is created with the ncdf4 package. Alternatively
## example data can be freely downloaded here: <https://wui.cmsaf.eu/>
library(ncdf4)
## create some (non-realistic) example data
lon <- seq(5, 15, 0.5)
lat <- seq(45, 55, 0.5)
time <- c(as.Date("2000-01-01"), as.Date("2001-02-01"))
origin <- as.Date("1983-01-01 00:00:00")
time <- as.numeric(difftime(time, origin, units = "hour"))
data1 <- array(250:350, dim = c(21, 21, 1))
data2 <- array(230:320, dim = c(21, 21, 1))
## create two simple example NetCDF files
x <- ncdim_def(name = "lon", units = "degrees_east", vals = lon)
y <- ncdim_def(name = "lat", units = "degrees_north", vals = lat)
t <- ncdim_def(name = "time", units = "hours since 1983-01-01 00:00:00",
vals = time[1], unlim = TRUE)
var1 <- ncvar_def("SIS", "W m-2", list(x, y, t), -1, prec = "short")
vars <- list(var1)
ncnew <- nc_create(file.path(tempdir(),"CMSAF_example_file_1.nc"), vars)
ncvar_put(ncnew, var1, data1)
ncatt_put(ncnew, "lon", "standard_name", "longitude", prec = "text")
ncatt_put(ncnew, "lat", "standard_name", "latitude", prec = "text")
nc_close(ncnew)
t <- ncdim_def(name = "time", units = "hours since 1983-01-01 00:00:00",
vals = time[2], unlim = TRUE)
ncnew <- nc_create(file.path(tempdir(),"CMSAF_example_file_2.nc"), vars)
ncvar_put(ncnew, var1, data2)
ncatt_put(ncnew, "lon", "standard_name", "longitude", prec = "text")
ncatt_put(ncnew, "lat", "standard_name", "latitude", prec = "text")
nc_close(ncnew)
## Cut a region and merge both example CM SAF NetCDF files into one
## output file. Get path information of working directory with getwd()
## command.
wd <- getwd()
cmsaf.cat(var = "SIS", infiles = c(file.path(tempdir(),
"CMSAF_example_file_1.nc"), file.path(tempdir(),"CMSAF_example_file_2.nc")),
outfile = file.path(tempdir(),"CMSAF_example_file_cat.nc"))
unlink(c(file.path(tempdir(),"CMSAF_example_file_1.nc"),
file.path(tempdir(),"CMSAF_example_file_2.nc"),
file.path(tempdir(),"CMSAF_example_file_cat.nc")))
[Package cmsafops version 1.4.0 Index]