whittaker.R 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529
  1. #' Filter Vegetation Index with Modified Whittaker Approach
  2. #'
  3. #' @description
  4. #' Use a modified Whittaker filter function (see References) from package
  5. #' \strong{ptw} to filter a vegetation index (VI) time serie of satellite data.
  6. #'
  7. #' @param vi \code{Raster*} or \code{character} filenames, sorted VI. Use
  8. #' \code{\link{preStack}} functionality to ensure the right input.
  9. #' @param w \code{Raster*} or \code{character} filenames. In case of MODIS
  10. #' composite, the sorted 'VI_Quality' layers.
  11. #' @param t \code{Raster*} or \code{character} filenames. In case of MODIS
  12. #' composite, the sorted 'composite_day_of_the_year' layers. If missing, the
  13. #' date is determined using \code{timeInfo}.
  14. #' @param timeInfo Output from \code{\link{orgTime}}.
  15. #' @param lambda \code{character} or \code{integer}. Yearly lambda value passed
  16. #' to \code{\link{whit2}}. If set as \code{character} (i.e.,
  17. #' \code{lambda = "600"}), it is not adapted to the time serie length, but used
  18. #' as a fixed value (see Details). High values = stiff/rigid spline.
  19. #' @param nIter \code{integer}. Number of iteration for the upper envelope
  20. #' fitting.
  21. #' @param outputAs \code{character}, organisation of output files.
  22. #' \code{"single"} (default) means each date one \code{RasterLayer};
  23. #' \code{"yearly"} a \code{RasterBrick} for each year, and \code{"one"} one
  24. #' \code{RasterBrick} for the entire time series.
  25. #' @param collapse \code{logical}. Collapse input data of multiple years into
  26. #' one single year before filtering.
  27. #' @param prefixSuffix \code{character}, file naming. Names are composed
  28. #' dot-separated: \code{paste0(prefixSuffix[1], "YYYDDD", lambda,
  29. #' refixSuffix[2], ".defaultFileExtension")}.
  30. #' @param outDirPath \code{character}, output path. Defaults to the current
  31. #' working directory.
  32. #' @param outlierThreshold \code{numeric} in the same unit as \code{vi}, used
  33. #' for outlier removal (see Details).
  34. #' @param mergeDoyFun Especially when using \code{collapse = TRUE}, multiple
  35. #' measurements for one day can be present. Here you can choose how those values
  36. #' are merged to one single value: \code{"max"} uses the highest value,
  37. #' \code{"mean"} or \code{"weighted.mean"} use the \code{\link{mean}} if no
  38. #' weighting \code{"w"} is available, and \code{\link{weighted.mean}} if it is.
  39. #' @param ... Arguments passed to \code{\link{writeRaster}} (except for
  40. #' \code{filename}).
  41. #'
  42. #' @return
  43. #' A Whittaker-smoothened \code{RasterStack}.
  44. #'
  45. #' @details
  46. #' The argument \code{lambda} is passed to \code{MODIS:::miwhitatzb1}. You can
  47. #' set it as yearly \code{lambda}, which means that it doesn't matter how long
  48. #' the input time serie is because \code{lambda} is adapted to it with:
  49. #' \code{lambda*('length of input timeserie in days'/365)}. The input length can
  50. #' differ from the output because of the \code{pillow} argument in
  51. #' \code{orgTime}. But it can also be set as \code{character} (i.e.,
  52. #' \code{lambda = "1000"}). In this case, the adaption to the time series length
  53. #' is not performed.\cr
  54. #'
  55. #' @references
  56. #' Modified Whittaker smoother, according to Atzberger & Eilers 2011
  57. #' International Journal of Digital Earth 4(5):365-386.\cr
  58. #' Implementation in R: Agustin Lobo 2012
  59. #'
  60. #' @note
  61. #' Currently tested on MODIS and Landsat data. Using M*D13, it is also possible
  62. #' to use the 'composite_day_of_the_year' and the 'VI_Quality' layers. Note that
  63. #' this function is currently under heavy development.
  64. #'
  65. #' @seealso
  66. #' \code{\link{smooth.spline.raster}}, \code{\link{raster}}.
  67. #'
  68. #' @author
  69. #' Matteo Mattiuzzi and Agustin Lobo
  70. #'
  71. #' @examples
  72. #' \dontrun{
  73. #' # The following function will download bit more than 1 year of MOD13A1 (~180mB) and therefore
  74. #' # take while to execute! Data will be downloaded to options("MODIS_localArcPath") and processed
  75. #' # to 'paste0(options("MODIS_outDirPath"),"fullCapa")'
  76. #' # You need to extract: 'vegetation index', 'VI_Quality layer', and 'composite day of the year',
  77. #' # this is expressed by the argument 'SDSstring'
  78. #' runGdal(product="MOD13A2",begin="2004340",extent="irland",end="2006020", job="fullCapa",
  79. #' SDSstring="101000000010")
  80. #' path <- paste0(options("MODIS_outDirPath"),"fullCapa")
  81. #'
  82. #' # the only obligatory dataset is the vegetatino index
  83. #' # get the 'vi' data in the source directory:
  84. #' vi <- preStack(path=path, pattern="*_NDVI.tif$")
  85. #'
  86. #' # "orgTime" detects timing information of the input data and generates based on the arguments
  87. #' # the output date information.
  88. #' # For spline functions (in general) the beginning and the end of the time series
  89. #' # is always problematic. So there is the argument "pillow" (default 75 days) that adds
  90. #' # (if available) some more layers on the two endings.
  91. #' timeInfo <- orgTime(vi,nDays=16,begin="2005001",end="2005365",pillow=40)
  92. #'
  93. #' # now re-run "preStack" with two differences, 'files' (output of the first 'preStack' call)
  94. #' # and the 'timeInfo'
  95. #' # Here only the data needed for the filtering is extracted:
  96. #' vi <- preStack(files=vi,timeInfo=timeInfo)
  97. #'
  98. #' whittaker.raster(vi,timeInfo=timeInfo,lambda=5000)
  99. #'
  100. #' # if the files are M*D13 you can use also Quality layers and the composite day of the year:
  101. #' wt <- preStack(path=path, pattern="*_VI_Quality.tif$", timeInfo=timeInfo)
  102. #' # can also be already stacked:
  103. #' inT <- preStack(path=path, pattern="*_composite_day_of_the_year.tif$", timeInfo=timeInfo)
  104. #'
  105. #' whittaker.raster(vi=vi, wt=wt, inT=inT, timeInfo=timeInfo, lambda=5000, overwrite=TRUE)
  106. #' }
  107. #'
  108. #' @name whittaker.raster
  109. #' @export whittaker.raster
  110. whittaker.raster <- function(vi, w=NULL, t=NULL, timeInfo = orgTime(vi), lambda = 5000, nIter= 3, outputAs="single", collapse=FALSE, prefixSuffix=c("MCD","ndvi"), outDirPath=".", outlierThreshold=NULL, mergeDoyFun="max", ...)
  111. {
  112. # debug
  113. # w=wt; t=inT; groupYears=TRUE; lambda = 5000; nIter= 3; outDirPath = "./"; collapse=FALSE; opt <- MODIS:::combineOptions();removeOutlier=FALSE;mergeDoyFun="max";opts$outDirPath <- MODIS:::setPath(outDirPath)
  114. # w=NULL; t=NULL; groupYears=TRUE; lambda = 5000; nIter= 5; outDirPath = "./SUB/";collapse=TRUE; opts <- MODIS:::combineOptions();removeOutlier=FALSE;mergeDoyFun="max"
  115. # opt <- list(bitShift=2,bitMask=15,threshold=6)
  116. # opt <- list(bitShift=2,bitMask=15)
  117. # opt <- list()
  118. opts <- combineOptions(...)
  119. outDirPath <- setPath(outDirPath)
  120. bitShift <- opts$bitShift
  121. bitMask <- opts$bitMask
  122. threshold <- opts$threshold
  123. dataFormat <- opts$dataFormat
  124. rasterOut <- toupper(raster::writeFormats())
  125. if(toupper(dataFormat) %in% rasterOut[,"name"]) {
  126. dataFormat <- getExtension(dataFormat)
  127. } else {
  128. stop("Unknown or unsupported data format: '", dataFormat, "'. Please run
  129. raster::writeFormats() (column 'name') for supported file types.\n")
  130. }
  131. minDat <- ifelse(is.null(opts$minDat), 3, opts$minDat) # 3 is very small!
  132. if (collapse) {
  133. if(timeInfo$call$nDays == "asIn")
  134. stop("Argument nDays = 'asIn' (passed to orgTime()) is not allowed when using collapse = TRUE.\n")
  135. fitt <- seq(as.numeric(format(min(timeInfo$outputLayerDates),"%j")),as.numeric(format(max(timeInfo$outputLayerDates),"%j")),by=timeInfo$call$nDays) + timeInfo$call$pillow
  136. } else {
  137. fitt <- timeInfo$outSeq
  138. }
  139. inlam <- lambda
  140. ## fixed lambda
  141. if (is.character(lambda)) {
  142. cat("Using fixed 'lambda': ", lambda, ".\n", sep = "")
  143. nameL <- "fL"
  144. ## yearly lambda
  145. } else {
  146. if (collapse) {
  147. lambda <- lambda * ((365 + 2 * timeInfo$call$pillow) / 365)
  148. cat("Yearly 'lambda' is:", inlam, "\nNow changed with lambda*((365+2*pillow)/365) to:",lambda,"\n")
  149. } else {
  150. lambda <- lambda * ((max(timeInfo$inSeq) - min(timeInfo$inSeq) - 1) / 365)
  151. cat("Yearly 'lambda' is: ", inlam, ".\n",
  152. "Now changed to lambda * ('length of input period in days' / 365): ",
  153. lambda, ".\n", sep = "")
  154. }
  155. nameL <- "yL"
  156. }
  157. if (is.character(inlam))
  158. inlam <- as.numeric(inlam)
  159. inlam <- round(inlam) #from here on used only in outputfilename
  160. lambda <- as.numeric(lambda)
  161. if (!inherits(vi, "Raster"))
  162. vi <- raster::stack(vi, quick = TRUE)
  163. if(!inherits(w, "Raster") & !is.null(w))
  164. w <- raster::stack(w, quick = TRUE)
  165. if(!inherits(t, "Raster") & !is.null(t))
  166. t <- raster::stack(t, quick = TRUE)
  167. if (is.null(opts$datatype))
  168. opts$datatype <- raster::dataType(vi[[1]])
  169. if (length(grep("FLT", opts$datatype)) > 0) {
  170. doround <- FALSE
  171. } else {
  172. doround <- TRUE
  173. }
  174. if (is.null(opts$overwrite)) {
  175. opts$overwrite <- FALSE
  176. }
  177. outputAs <- tolower(outputAs)
  178. if (collapse) # use only DOY obmitt year (all years into one)
  179. {
  180. d <- sprintf("%03d",seq(as.numeric(format(min(timeInfo$outputLayerDates),"%j")),as.numeric(format(max(timeInfo$outputLayerDates),"%j")),by=timeInfo$call$nDays))
  181. if(outputAs=="single")
  182. {
  183. oname <- paste0(outDirPath,prefixSuffix[1],".",d,".",prefixSuffix[2],dataFormat)
  184. b <- vector(mode="list",length=length(oname))
  185. for(a in seq_along(oname))
  186. {
  187. b[[a]] <- raster(vi)
  188. b[[a]] <- writeStart(b[[a]], filename=oname[a], datatype=opts$datatype, overwrite=opts$overwrite)
  189. }
  190. } else
  191. {
  192. b <- list()
  193. b[[1]] <- writeStart(brick(raster(vi),nl=length(d)), filename=oname, datatype=opts$datatype, overwrite=opts$overwrite)
  194. names(b[[1]]) <- paste0("doy",d)
  195. }
  196. } else if (outputAs=="yearly")
  197. {
  198. b <- vector(mode="list",length=length(unique(format(timeInfo$outputLayerDates,"%Y"))))
  199. for (a in seq_along(unique(format(timeInfo$outputLayerDates,"%Y"))))
  200. {
  201. y <- unique(format(timeInfo$outputLayerDates,"%Y"))[a]
  202. oname <- paste0(outDirPath,prefixSuffix[1],".year",y,".",nameL,inlam,".",prefixSuffix[2],dataFormat)
  203. names <- timeInfo$outputLayerDates[format(timeInfo$outputLayerDates,"%Y")==y]
  204. b[[a]] <- brick(raster(vi),nl=as.integer(sum(format(timeInfo$outputLayerDates,"%Y")==y)), values=FALSE)
  205. b[[a]] <- writeStart(b[[a]], filename=oname, datatype=opts$datatype, overwrite=opts$overwrite)
  206. names(b[[a]]) <- names
  207. }
  208. } else if(outputAs=="one")
  209. {
  210. y <- unique(format(timeInfo$outputLayerDates,"%Y"))
  211. oname <- paste0(outDirPath,prefixSuffix[1],"_from",paste0(y,collapse="to"),".",nameL,inlam,".",prefixSuffix[2],dataFormat)
  212. b <- list()
  213. b[[1]] <- brick(raster(vi),nl=length(fitt), values=FALSE)
  214. b[[1]] <- writeStart(b[[1]], filename=oname, datatype=opts$datatype, overwrite=opts$overwrite)
  215. names(b[[a]]) <- timeInfo$outputLayerDates
  216. } else if (outputAs=="single")
  217. {
  218. d <- sort(format(timeInfo$outputLayerDates,"%Y%j"))
  219. oname <- paste0(outDirPath,prefixSuffix[1],".",d,".",nameL,inlam,".",prefixSuffix[2],dataFormat)
  220. b <- vector(mode="list",length=length(oname))
  221. for(a in seq_along(oname))
  222. {
  223. b[[a]] <- raster::raster(vi)
  224. b[[a]] <- writeStart(b[[a]], filename=oname[a], datatype=opts$datatype, overwrite=opts$overwrite)
  225. }
  226. }
  227. tr <- raster::blockSize(vi)
  228. cat("Data is in, start processing!\n")
  229. if (mergeDoyFun == "max") {
  230. mergeFun <- unifyDoubleMX
  231. } else if (mergeDoyFun == "weighted.mean" | mergeDoyFun == "mean") {
  232. mergeFun <- unifyDoubleWM
  233. }
  234. clFun <- function(l) {
  235. val <- raster::getValues(vi, row = tr$row[l], nrows = tr$nrows[l])
  236. val <- t(val)
  237. mtrdim <- dim(val)
  238. set0 <- matrix(FALSE,nrow = mtrdim[1], ncol = mtrdim[2])
  239. set0[is.na(val)] <- TRUE
  240. ## if 'VI_Quality' is supplied:
  241. if (!is.null(w)) {
  242. wtu <- raster::getValues(w, row = tr$row[l], nrows = tr$nrows[l])
  243. # is it not a weight info [0-1]?
  244. if (max(wtu, na.rm = TRUE) > 1) {
  245. if(is.null(bitShift) | is.null(bitMask)) {
  246. # try to detect VI usefulness layer
  247. bits <- detectBitInfo(vi, "VI usefulness", warn = FALSE)
  248. bitShift <- bits$bitShift
  249. bitMask <- bits$bitMask
  250. }
  251. if(is.null(bitShift) | is.null(bitMask))
  252. stop("Could not extract 'bits' for weighting from this product. ",
  253. "Use '?makeWeights' function to generate weights manually!")
  254. wtu <- makeWeights(wtu, bitShift = bitShift, bitMask = bitMask, threshold = threshold, decodeOnly = FALSE)
  255. }
  256. wtu <- t(wtu)
  257. set0[wtu==0] <- TRUE
  258. ## else if 'VI_Quality' is not supplied, then weight = 1:
  259. } else {
  260. wtu <- matrix(1, nrow = mtrdim[1], ncol = mtrdim[2])
  261. }
  262. if (inherits(t, "Raster")) {
  263. inTu <- raster::getValues(t, row = tr$row[l], nrows = tr$nrows[l])
  264. inTu <- t(inTu)
  265. set0[is.na(inTu)] <- TRUE
  266. set0[inTu <= 0] <- TRUE
  267. t0 <- min(timeInfo$inDoys[1]) - 1
  268. if (!collapse) {
  269. inTu <- t(repDoy(t(inTu), layerDate = timeInfo, bias = -t0))
  270. }
  271. inTu[set0] <- 0
  272. } else {
  273. if (collapse) {
  274. inTu <- matrix(timeInfo$inDoys,nrow=length(timeInfo$inDoys),ncol=mtrdim[2])
  275. } else {
  276. inTu <- matrix(timeInfo$inSeq,nrow=length(timeInfo$inSeq),ncol=mtrdim[2])
  277. }
  278. }
  279. # the entire info to use or not a pix is in "wtu"
  280. wtu[set0] <- 0
  281. val[set0] <- 0
  282. out <- matrix(NA, nrow = length(fitt), ncol = mtrdim[2])
  283. if (!is.null(outlierThreshold)) {
  284. kickOutlier <- function(vals, weights, lambda, threshold) {
  285. fTS <- ptw::whit2(vals, w = weights, lambda = lambda)
  286. weights[weights==1][abs(vals[weights==1]-fTS[weights==1]) > threshold] <- 0
  287. return(weights)
  288. }
  289. } else {
  290. # if is.null(outlierThreshold) generate a fake function to avoid a per pixel "if"
  291. kickOutlier <- function(vals, weights, lambda, threshold) {
  292. return(weights)
  293. }
  294. }
  295. if (collapse)
  296. {
  297. vec0 <- rep(0,365 + (2*timeInfo$call$pillow) + 30) # add a save length of data (because layer doy + effectice composite doy)
  298. } else
  299. {
  300. vec0 <- rep(0,max(timeInfo$inSeq,timeInfo$outSeq) - min(timeInfo$inSeq,timeInfo$outSeq) - 1 + 30)
  301. }
  302. # minimum "minDat" input values for filtering
  303. Cvec <- (colSums(wtu > 0) >= minDat)
  304. Cvec <- (1:mtrdim[2])[Cvec]
  305. ind <- inTu > 0
  306. win <- options("warn")
  307. options(warn=-1)
  308. for (u in Cvec)
  309. {
  310. index <- ind[,u]
  311. use <- mergeFun(vx=val[index,u],wx=wtu[index,u],tx=inTu[index,u])
  312. valVec <- wtVec <- vec0
  313. if(!collapse)
  314. {
  315. valVec[use$tx] <- use$vx
  316. wtVec[use$tx] <- use$wx
  317. } else
  318. {
  319. newOrder <- doCollapse(tx=use$tx,pillow=timeInfo$call$pillow)
  320. valVec[newOrder$sequence] <- use$vx[newOrder$order]
  321. wtVec[newOrder$sequence] <- use$wx[newOrder$order]
  322. }
  323. wtVec <- kickOutlier(vals=valVec,weights=wtVec,lambda=lambda,threshold=outlierThreshold)
  324. #plot(valVec,ylim=c(-1000,9000))
  325. for(i in 1:nIter)
  326. {
  327. fTS <- ptw::whit2(valVec,w=wtVec,lambda=lambda)
  328. valVec[valVec < fTS] <- fTS[valVec < fTS]
  329. }
  330. out[,u] <- fTS[fitt]
  331. #lines(fTS,col=2)
  332. }
  333. options(warn=win$warn)
  334. out[,colSums(abs(out))==0] <- NA
  335. return(t(out))
  336. }
  337. for (i in seq_along(tr$row)) {
  338. res <- clFun(i)
  339. if (doround)
  340. res <- round(res)
  341. b <- writeValuesMODIS(b, res, tr$row[i], timeInfo, collapse, outputAs)
  342. }
  343. writeStopMODIS(b,timeInfo,outputAs,collapse)
  344. return(raster::stack(b))
  345. }
  346. unifyDoubleWM <- function(vx,wx,tx)
  347. {
  348. tx <- as.numeric(tx)
  349. vx <- as.numeric(vx)
  350. wx <- as.numeric(wx)
  351. double <- tx[duplicated(tx)]
  352. if(length(double)>0)
  353. {
  354. double <- unique(double)
  355. for(i in seq_along(double))
  356. {
  357. inx <- which(tx==double[i])
  358. vx[inx[1]] <- weighted.mean(vx[inx],w=wx[inx])
  359. wx[inx[1]] <- max(wx[inx])
  360. vx <- vx[-inx[-1]]
  361. wx <- wx[-inx[-1]]
  362. tx <- tx[-inx[-1]]
  363. }
  364. }
  365. list(vx=vx,wx=wx,tx=tx)
  366. }
  367. unifyDoubleMX <- function(vx,wx,tx)
  368. {
  369. tx <- as.numeric(tx)
  370. vx <- as.numeric(vx)
  371. wx <- as.numeric(wx)
  372. double <- tx[duplicated(tx)]
  373. if(length(double)>0)
  374. {
  375. double <- unique(double)
  376. for(i in seq_along(double))
  377. {
  378. inx <- which(tx==double[i])
  379. mx <- which.max(wx[inx])
  380. vx <- vx[-inx[-mx]]
  381. wx <- wx[-inx[-mx]]
  382. tx <- tx[-inx[-mx]]
  383. }
  384. }
  385. list(vx=vx,wx=wx,tx=tx)
  386. }
  387. doCollapse <- function(tx,pillow)
  388. {
  389. ord <- order(tx)
  390. txS <- tx[ord]
  391. t0 <- 365 - pillow
  392. tS <- ord[txS >= t0]
  393. tE <- ord[txS <= pillow]
  394. s0 <- txS[txS >= t0] - t0
  395. s1 <- txS + pillow
  396. s2 <- txS[txS <= pillow] + 365 + pillow
  397. list(order=c(tS,ord,tE),sequence=c(s0,s1,s2)+1)
  398. }
  399. #####################################################
  400. writeValuesMODIS <- function(b,val,row,timeInfo,collapse,outputAs)
  401. {
  402. if(collapse)
  403. {
  404. d <- seq_along(sprintf("%03d",seq(as.numeric(format(min(timeInfo$outputLayerDates),"%j")),as.numeric(format(max(timeInfo$outputLayerDates),"%j")),by=timeInfo$call$nDays)))
  405. } else
  406. {
  407. d <- seq_along(b)
  408. }
  409. if(outputAs=="single")
  410. {
  411. for (a in seq_along(d))
  412. {
  413. b[[a]] <- writeValues(b[[a]], val[,a], row)
  414. }
  415. } else if(outputAs=="one")
  416. {
  417. b[[1]] <- writeValues(b[[1]], val, row)
  418. } else
  419. {
  420. for (a in seq_along(d))
  421. {
  422. y <- unique(format(timeInfo$outputLayerDates,"%Y"))[a]
  423. b[[a]] <- writeValues(b[[a]], val[,format(timeInfo$outputLayerDates,"%Y")==y], row)
  424. }
  425. }
  426. return(b)
  427. }
  428. writeStopMODIS <- function(b,timeInfo,outputAs,collapse)
  429. {
  430. for (a in seq_along(b))
  431. {
  432. b[[a]] <- writeStop(b[[a]])
  433. nam <- filename(b[[a]])
  434. extension(nam) <- ""
  435. if (collapse & outputAs!="single")
  436. {
  437. write.table(x=unique(format(timeInfo$outputLayerDates,"%j")),
  438. file=nam, col.names=FALSE, row.names=FALSE)
  439. } else if(outputAs=="one")
  440. {
  441. write.table(x=timeInfo$outputLayerDates, file=nam,
  442. col.names=FALSE, row.names=FALSE)
  443. } else if (outputAs=="yearly")
  444. {
  445. y <- unique(format(timeInfo$outputLayerDates,"%Y"))
  446. for (ax in seq_along(y))
  447. {
  448. ind <- format(timeInfo$outputLayerDates,"%Y")==y[ax]
  449. write.table(x=timeInfo$outputLayerDates[ind], file=nam,
  450. col.names=FALSE, row.names=FALSE)
  451. }
  452. }
  453. }
  454. }