Skip to content

Commit

Permalink
return sftime object; addresses #104
Browse files Browse the repository at this point in the history
  • Loading branch information
edzer committed Jun 18, 2022
1 parent 3bf5981 commit d5ba538
Show file tree
Hide file tree
Showing 6 changed files with 72 additions and 11 deletions.
6 changes: 2 additions & 4 deletions DESCRIPTION
Original file line number Diff line number Diff line change
Expand Up @@ -10,10 +10,8 @@ Authors@R: c(person(given = "Edzer",
person("Benedikt", "Graeler", role = "aut"))
Description: Variogram modelling; simple, ordinary and universal point or block (co)kriging; spatio-temporal kriging; sequential Gaussian or indicator (co)simulation; variogram and variogram map plotting utility functions; supports sf and stars.
Depends: R (>= 2.10)
Imports: utils, stats, graphics, methods, lattice, sp (>= 0.9-72), zoo,
spacetime (>= 1.2-7), FNN
Suggests: fields, maps, mapdata, maptools, rgdal (>= 0.5.2), rgeos, sf
(>= 0.7-2), sftime, stars (>= 0.3-0), xts, raster, future, future.apply
Imports: utils, stats, graphics, methods, lattice, sp (>= 0.9-72), zoo, sf (>= 0.7-2), sftime, spacetime (>= 1.2-8), stars, FNN
Suggests: fields, maps, mapdata, maptools, rgdal (>= 0.5.2), rgeos, xts, raster, future, future.apply
License: GPL (>= 2.0)
URL: https://github.com/r-spatial/gstat/
Encoding: UTF-8
Expand Down
2 changes: 1 addition & 1 deletion NEWS.md
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
# version 2.1-0

* import `sftime`; modify `krigeST()` variogram functions to accept `sftime` objects as alternative to `STI` or `STIDF`
* import `sftime`; modify `krigeST()` variogram functions to accept `sftime` objects for `data` (as alternative to `STI` or `STIDF`), and `stars` or `sftime` objects for `newdata`; #108 with great help from @henningte

* fix Gaussian cosimulation, probably introduced in 2016, #106

Expand Down
16 changes: 12 additions & 4 deletions R/krigeST.R
Original file line number Diff line number Diff line change
Expand Up @@ -54,6 +54,7 @@ krigeST <- function(formula, data, newdata, modelList, beta, y, ...,
computeVar = FALSE, fullCovariance = FALSE,
bufferNmax=2, progress=TRUE) {
stopifnot(inherits(modelList, "StVariogramModel") || is.function(modelList))
to_sftime = FALSE
return_stars = if (inherits(data, c("stars"))) {
if (!requireNamespace("sf", quietly = TRUE))
stop("sf required: install that first") # nocov
Expand All @@ -73,8 +74,10 @@ krigeST <- function(formula, data, newdata, modelList, beta, y, ...,
newdata$._dummy = 0.
newdata = as(newdata, "STFDF")
}
if (inherits(newdata, "sftime"))
if (inherits(newdata, "sftime")) {
to_sftime = TRUE
newdata = as(newdata, "STIDF")
}
TRUE
} else {
if (!identical(data@sp@proj4string@projargs, newdata@sp@proj4string@projargs))
Expand All @@ -83,7 +86,6 @@ krigeST <- function(formula, data, newdata, modelList, beta, y, ...,
}
stopifnot(inherits(data, c("STF", "STS", "STI", "sftime")) &&
inherits(newdata, c("STF", "STS", "STI", "sftime")))
to_sftime = inherits(data, "sftime") # deal with later
if (inherits(data, "sftime"))
data = as(data, "STIDF")
if (inherits(newdata, "sftime"))
Expand Down Expand Up @@ -113,7 +115,10 @@ krigeST <- function(formula, data, newdata, modelList, beta, y, ...,

if (return_stars) {
ret$._dummy = NULL
stars::st_as_stars(as(ret, "STFDF"))
if (to_sftime)
sftime::st_as_sftime(ret)
else
stars::st_as_stars(as(ret, "STFDF"))
} else
ret
} else {
Expand All @@ -128,8 +133,11 @@ krigeST <- function(formula, data, newdata, modelList, beta, y, ...,
if (!fullCovariance) {
ret = addAttrToGeom(geometry(newdata), df)
if (return_stars) {
ret = stars::st_as_stars(as(ret, "STFDF"))
ret$._dummy = NULL
if (to_sftime)
ret = sftime::st_as_sftime(ret)
else
ret = stars::st_as_stars(as(ret, "STFDF"))
}
ret
} else
Expand Down
6 changes: 6 additions & 0 deletions demo/sftime.R
Original file line number Diff line number Diff line change
Expand Up @@ -53,3 +53,9 @@ stplot(locKrig[,,"var1.pred"], col.regions=bpy.colors(), scales=list(draw=T))
plot(locKrig_sft[1], col = sf.colors(), breaks = "equal")
stplot(locKrig[,,"var1.var"], col.regions=bpy.colors(), scales=list(draw=T))
plot(locKrig_sft[2], col = sf.colors(), breaks = "equal")

st$foo = 0
st_as_sf(st, long = TRUE) |> st_as_sftime() -> st.sftime
locKrig_sft <- krigeST(z~1, sft, st.sftime, sumMetricModel, nmax=20, computeVar = T)
plot(locKrig_sft["var1.pred"])

43 changes: 43 additions & 0 deletions inst/ChangeLog
Original file line number Diff line number Diff line change
@@ -1,3 +1,46 @@
2022-06-12 Edzer Pebesma <[email protected]>

* R/krigeST.R: remove ._dummy variable; #108

2022-06-05 Edzer Pebesma <[email protected]>

* demo/00Index: update index

2022-06-05 Edzer Pebesma <[email protected]>

* R/krigeST.R, demo/sftime.R: also handle stars with 0 arrays as
newdata

2022-06-02 Edzer Pebesma <[email protected]>

* demo/sftime.R: add demo script

2022-06-02 Edzer Pebesma <[email protected]>

* DESCRIPTION, R/krigeST.R, R/stVariogramModels.R, R/variogramST.R,
inst/ChangeLog, tests/stars.Rout.save: modifications to use sftime &
stars.

2022-05-31 Edzer Pebesma <[email protected]>

* .Rbuildignore: tidy

2022-05-31 Edzer Pebesma <[email protected]>

* NEWS.md: add news

2022-05-31 Edzer Pebesma <[email protected]>

* src/mtrx.c: tidy

2022-05-01 Edzer Pebesma <[email protected]>

* .github/workflows/tic.yml: GA, see #107

2022-04-29 Edzer Pebesma <[email protected]>

* src/sim.c: fixes #106

2022-04-12 Edzer Pebesma <[email protected]>

* R/stVariogramModels.R: replace | with || in if() calls
Expand Down
10 changes: 8 additions & 2 deletions man/krigeST.Rd
Original file line number Diff line number Diff line change
Expand Up @@ -64,6 +64,11 @@ krigeSTTg(formula, data, newdata, modelList, y, nmax=Inf, stAni=NULL,
particular methods for ordinary spatio-temporal (ST) kriging. In particular,
it does not support block kriging or kriging in a distance-based
neighbourhood, and does not provide simulation.

If \code{data} is of class \code{sftime}, then \code{newdata} MUST be
of class \code{stars} or \code{sftime}, i.e. mixing form old-style
classes (package spacetime) and new-style classes (sf, stars, sftime)
is not supported.
}

\value{
Expand All @@ -73,8 +78,9 @@ krigeSTTg(formula, data, newdata, modelList, y, nmax=Inf, stAni=NULL,
}

\references{
Spatio-Temporal Geostatistics using gstat. Benedikt Graeler,
Edzer Pebesma, Gerard Heuvelink. The R Journal, accepted.
Benedikt Graeler, Edzer Pebesma, Gerard Heuvelink. Spatio-Temporal
Geostatistics using gstat. The R Journal 8(1), 204--218.
\url{https://journal.r-project.org/archive/2016/RJ-2016-014/index.html}

N.A.C. Cressie, 1993, Statistics for Spatial Data,
Wiley.
Expand Down

0 comments on commit d5ba538

Please sign in to comment.