Title: | Tracking Data |
---|---|
Description: | Access and manipulate spatial tracking data, with straightforward coercion from and to other formats. Filter for speed and create time spent maps from tracking data. There are coercion methods to convert between 'trip' and 'ltraj' from 'adehabitatLT', and between 'trip' and 'psp' and 'ppp' from 'spatstat'. Trip objects can be created from raw or grouped data frames, and from types in the 'sp', sf', 'amt', 'trackeR', 'mousetrap', and other packages, Sumner, MD (2011) <https://figshare.utas.edu.au/articles/thesis/The_tag_location_problem/23209538>. |
Authors: | Michael D. Sumner [aut, cre], Sebastian Luque [ctb], Anthony Fischbach [ctb], Tomislav Hengl [ctb] |
Maintainer: | Michael D. Sumner <[email protected]> |
License: | GPL-3 |
Version: | 1.10.0.9010 |
Built: | 2024-11-03 05:29:24 UTC |
Source: | https://github.com/trackage/trip |
Functions for accessing and manipulating spatial data for animal tracking, with straightforward coercion from and to other formats. Filter for speed and create time spent maps from animal track data. There are coercion methods to convert between 'trip' and 'ltraj' from 'adehabitatLT', and between 'trip' and 'psp' and 'ppp' from 'spatstat'. Trip objects can be created from raw or grouped data frames, and from types in the 'sp', 'sf', 'amt', 'trackeR', and other packages.
Maintainer: Michael D. Sumner [email protected]
Other contributors:
Sebastian Luque [contributor]
Anthony Fischbach [contributor]
Tomislav Hengl [contributor]
Useful links:
Duplicated DateTime values within ID are adjusted forward (recursively) by one second until no duplicates are present. This is considered reasonable way of avoiding the nonsensical problem of duplicate times.
adjust.duplicateTimes(time, id)
adjust.duplicateTimes(time, id)
time |
vector of DateTime values |
id |
vector of ID values, matching DateTimes that are assumed sorted within ID |
This function is used to remove duplicate time records in animal track data, rather than removing the record completely.
The adjusted DateTime vector is returned.
I have no idea what goes on at CLS when they output data that are either not ordered by time or have duplicates. If this problem exists in your data it's probably worth finding out why.
## DateTimes with a duplicate within ID tms <- Sys.time() + c(1:6, 6, 7:10) *10 id <- rep("a", length(tms)) range(diff(tms)) ## duplicate record is now moved one second forward tms.adj <- adjust.duplicateTimes(tms, id) range(diff(tms.adj))
## DateTimes with a duplicate within ID tms <- Sys.time() + c(1:6, 6, 7:10) *10 id <- rep("a", length(tms)) range(diff(tms)) ## duplicate record is now moved one second forward tms.adj <- adjust.duplicateTimes(tms, id) range(diff(tms.adj))
Assign numeric values for Argos "class" by matching the levels available to given numbers. An adjustment is made to allow sigma to be specified in kilometres, and the values returned are the approximate values for longlat degrees. It is assumed that the levels are part of an "ordered" factor from least precise to most precise.
argos.sigma(x, sigma = c(100, 80, 50, 20, 10, 4, 2), adjust = 111.12)
argos.sigma(x, sigma = c(100, 80, 50, 20, 10, 4, 2), adjust = 111.12)
x |
factor of Argos location quality "classes" |
sigma |
numeric values (by default in kilometres) |
adjust |
a numeric adjustment to convert from kms to degrees |
The available levels in Argos are levels=c("Z", "B", "A", "0", "1",
"2", "3")
.
The actual sigma values given by default are (as far as can be determined) a reasonable stab at what Argos believes.
Numeric values for given levels.
cls <- ordered(sample(c("Z", "B", "A", "0", "1", "2", "3"), 30, replace=TRUE), levels=c("Z", "B", "A", "0", "1", "2", "3")) argos.sigma(cls)
cls <- ordered(sample(c("Z", "B", "A", "0", "1", "2", "3"), 30, replace=TRUE), levels=c("Z", "B", "A", "0", "1", "2", "3")) argos.sigma(cls)
Coercing trip
objects to other classes.
Function to create a SpatialLinesDataFrame from a trip object, resulting in a line segment for each implicit segment along the tracks. The object stores the start and end times, duration and the ID of the segment.
## S3 method for class 'trip' as.ppp(X, ..., fatal) ## S3 method for class 'trip' as.psp(x, ..., from, to) as.track_xyt.trip(x, ..., from, to) explode(x, ...)
## S3 method for class 'trip' as.ppp(X, ..., fatal) ## S3 method for class 'trip' as.psp(x, ..., from, to) as.track_xyt.trip(x, ..., from, to) explode(x, ...)
X |
|
... |
reserved for future methods |
fatal |
Logical value, see Details of |
x |
|
from |
see |
to |
See |
ppp object
psp object
SpatialLinesDataFrame
SpatialLinesDataFrame object with each individual line segment identified by start/end time and trip ID
d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) sp::coordinates(d) <- ~x+y tr <- trip(d, c("tms", "id")) as(tr, "ppp") d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d, c("tms", "id")) as(tr, "psp") as.psp(tr) d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d) spldf <- explode(tr) summary(tr)
d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) sp::coordinates(d) <- ~x+y tr <- trip(d, c("tms", "id")) as(tr, "ppp") d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d, c("tms", "id")) as(tr, "psp") as.psp(tr) d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d) spldf <- explode(tr) summary(tr)
trip
objectsCoercing objects to trip
class
as.trip(x, ...)
as.trip(x, ...)
x |
ltr ltraj object |
... |
Arguments passed to other methods. Ignored for |
S4 trip object
signature(from="ltraj", to="trip")
signature(x="ltraj")
d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d) if (require(adehabitatLT)) { l <- as(tr, "ltraj") ltraj2trip(l) as.trip(l) }
d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d) if (require(adehabitatLT)) { l <- as(tr, "ltraj") ltraj2trip(l) as.trip(l) }
Split trip events within a single object into exact time boundaries, adding interpolated coordinates as required.
## S3 method for class 'trip' cut(x, breaks, ...)
## S3 method for class 'trip' cut(x, breaks, ...)
x |
A trip object. |
breaks |
A character string such as the |
... |
Unused arguments. |
Motion between boundaries is assumed linear and extra coordinates are added at the cut points.
This function was completely rewritten in version 1.1-20.
list of S4 trip objects, each with aligned boundaries in time based on cutting the input into intervals
A list of trip objects, named by the time boundary in which they lie.
Michael D. Sumner and Sebastian Luque
See also tripGrid
.
set.seed(66) d <- data.frame(x=1:100, y=rnorm(100, 1, 10), tms= as.POSIXct(as.character(Sys.time()), tz = "GMT") + c(seq(10, 1000, length=50), seq(100, 1500, length=50)), id=gl(2, 50)) sp::coordinates(d) <- ~x+y tr <- trip(d, c("tms", "id")) cut(tr, "200 sec") bound.dates <- seq(min(tr$tms) - 1, max(tr$tms) + 1, length=5) trip.list <- cut(tr, bound.dates) bb <- sp::bbox(tr) cn <- c(20, 8) g <- sp::GridTopology(bb[, 1], apply(bb, 1, diff) / (cn - 1), cn) tg <- tripGrid(tr, grid=g) tg <- sp::as.image.SpatialGridDataFrame(tg) tg$x <- tg$x - diff(tg$x[1:2]) / 2 tg$y <- tg$y - diff(tg$y[1:2]) / 2 op <- par(mfcol=c(4, 1)) for (i in 1:length(trip.list)) { plot(sp::coordinates(tr), pch=16, cex=0.7) title(names(trip.list)[i], cex.main=0.9) lines(trip.list[[i]]) abline(h=tg$y, v=tg$x, col="grey") image(tripGrid(trip.list[[i]], grid=g), interpolate=FALSE, col=c("white", grey(seq(0.2, 0.7, length=256))),add=TRUE) abline(h=tg$y, v=tg$x, col="grey") lines(trip.list[[i]]) points(trip.list[[i]], pch=16, cex=0.7) } par(op) print("you may need to resize the window to see the grid data") cn <- c(200, 80) g <- sp::GridTopology(bb[, 1], apply(bb, 1, diff) / (cn - 1), cn) tg <- tripGrid(tr, grid=g) tg <- sp::as.image.SpatialGridDataFrame(tg) tg$x <- tg$x - diff(tg$x[1:2]) / 2 tg$y <- tg$y - diff(tg$y[1:2]) / 2 op <- par(mfcol=c(4, 1)) for (i in 1:length(trip.list)) { plot(sp::coordinates(tr), pch=16, cex=0.7) title(names(trip.list)[i], cex.main=0.9) image(tripGrid(trip.list[[i]], grid=g, method="density", sigma=1), interpolate=FALSE, col=c("white", grey(seq(0.2, 0.7, length=256))), add=TRUE) lines(trip.list[[i]]) points(trip.list[[i]], pch=16, cex=0.7) } par(op) print("you may need to resize the window to see the grid data") data("walrus818", package = "trip") library(lubridate) walrus_list <- cut(walrus818, seq(floor_date(min(walrus818$DataDT), "month"), ceiling_date(max(walrus818$DataDT), "month"), by = "1 month")) g <- rasterize(walrus818) * NA_real_ stk <- raster::stack(lapply(walrus_list, rasterize, grid = g)) st <- raster::aggregate(stk, fact = 4, fun = sum, na.rm = TRUE) st[!st > 0] <- NA_real_ plot(st, col = oc.colors(52))
set.seed(66) d <- data.frame(x=1:100, y=rnorm(100, 1, 10), tms= as.POSIXct(as.character(Sys.time()), tz = "GMT") + c(seq(10, 1000, length=50), seq(100, 1500, length=50)), id=gl(2, 50)) sp::coordinates(d) <- ~x+y tr <- trip(d, c("tms", "id")) cut(tr, "200 sec") bound.dates <- seq(min(tr$tms) - 1, max(tr$tms) + 1, length=5) trip.list <- cut(tr, bound.dates) bb <- sp::bbox(tr) cn <- c(20, 8) g <- sp::GridTopology(bb[, 1], apply(bb, 1, diff) / (cn - 1), cn) tg <- tripGrid(tr, grid=g) tg <- sp::as.image.SpatialGridDataFrame(tg) tg$x <- tg$x - diff(tg$x[1:2]) / 2 tg$y <- tg$y - diff(tg$y[1:2]) / 2 op <- par(mfcol=c(4, 1)) for (i in 1:length(trip.list)) { plot(sp::coordinates(tr), pch=16, cex=0.7) title(names(trip.list)[i], cex.main=0.9) lines(trip.list[[i]]) abline(h=tg$y, v=tg$x, col="grey") image(tripGrid(trip.list[[i]], grid=g), interpolate=FALSE, col=c("white", grey(seq(0.2, 0.7, length=256))),add=TRUE) abline(h=tg$y, v=tg$x, col="grey") lines(trip.list[[i]]) points(trip.list[[i]], pch=16, cex=0.7) } par(op) print("you may need to resize the window to see the grid data") cn <- c(200, 80) g <- sp::GridTopology(bb[, 1], apply(bb, 1, diff) / (cn - 1), cn) tg <- tripGrid(tr, grid=g) tg <- sp::as.image.SpatialGridDataFrame(tg) tg$x <- tg$x - diff(tg$x[1:2]) / 2 tg$y <- tg$y - diff(tg$y[1:2]) / 2 op <- par(mfcol=c(4, 1)) for (i in 1:length(trip.list)) { plot(sp::coordinates(tr), pch=16, cex=0.7) title(names(trip.list)[i], cex.main=0.9) image(tripGrid(trip.list[[i]], grid=g, method="density", sigma=1), interpolate=FALSE, col=c("white", grey(seq(0.2, 0.7, length=256))), add=TRUE) lines(trip.list[[i]]) points(trip.list[[i]], pch=16, cex=0.7) } par(op) print("you may need to resize the window to see the grid data") data("walrus818", package = "trip") library(lubridate) walrus_list <- cut(walrus818, seq(floor_date(min(walrus818$DataDT), "month"), ceiling_date(max(walrus818$DataDT), "month"), by = "1 month")) g <- rasterize(walrus818) * NA_real_ stk <- raster::stack(lapply(walrus_list, rasterize, grid = g)) st <- raster::aggregate(stk, fact = 4, fun = sum, na.rm = TRUE) st[!st > 0] <- NA_real_ plot(st, col = oc.colors(52))
A convenience function, that removes duplicate rows, sorts by the date-times within ID, and removes duplicates from a data frame or SpatialPointsDataFrame.
forceCompliance(x, tor)
forceCompliance(x, tor)
x |
|
tor |
character vector of names of date-times and trip ID columns |
data.frame
or
SpatialPointsDataFrame-class
.
It's really important that data used are of a given quality, but this function makes the most common trip problems easy to apply.
This function returns a distance from a given 'home' coordinate for each individual trip.
Use the home
argument to provide a single, common 2-element (x,y or lon,lat) coordinate. If home
is NULL
(the default), then each individual trip's first location is used.
homedist(x, home = NULL)
homedist(x, home = NULL)
x |
trip object |
home |
see details |
numeric vector of distances in km (for longlat), or in the units of the trip's projection
Calculate great circle intermediate points on longitude, latitude input vectors. A spherical model is used, from the geosphere package.
interp_equal(x, distance = NULL, duration = NULL)
interp_equal(x, distance = NULL, duration = NULL)
x |
trip object |
distance |
optional minimum distance (metres) between interpolated points |
duration |
optional minimum duration (seconds) between interpolated points |
For the result to be sensible, the input must either be in longitude/latitude, or be in a projection and have a valid CRS. Great circle movement is assumed, there's no way to use this to interpolate equal-distance in the native projection.
If no input distance
or duration
is provided a default is used of 15 points between each input point.
if both distance
AND duration
is provided, distance
is ignored.
Note, the original implementation of this function was called 'interpequal()', and was used for time spent calculations. The functionality is now provided by the traipse package.
S4 trip object with interpolated new locations based on distance or duration parameters
Sensible defaults are assumed, to match the extents of data to a manageable grid.
makeGridTopology( obj, cells.dim = c(100, 100), xlim = NULL, ylim = NULL, buffer = 0, cellsize = NULL, adjust2longlat = FALSE )
makeGridTopology( obj, cells.dim = c(100, 100), xlim = NULL, ylim = NULL, buffer = 0, cellsize = NULL, adjust2longlat = FALSE )
obj |
any Spatial object, or other object for which |
cells.dim |
the number of cells of the grid, x then y |
xlim |
x limits of the grid |
ylim |
y limits of the grid |
buffer |
proportional size of the buffer to add to the grid limits |
cellsize |
pixel cell size |
adjust2longlat |
assume cell size is in kilometres and provide simple adjustment for earth-radius cells at the north-south centre of the grid |
Approximations for kilometres in longlat can be made using cellsize
and adjust2longlat
.
S4 class GridTopology with properties set variously from input parameters
Generate ocean colour colours, using the SeaWiFS scheme
oc.theme(x = 50) oc.colors(n)
oc.theme(x = 50) oc.colors(n)
x |
Number of colours to generate as part of a theme |
n |
Number of colours to generate |
This is a high-contrast palette, log-scaled originally for ocean chlorophyll.
A set of colours or a theme object.
Similar functions in sp spplot
,
bpy.colors
oc.colors(10) library(lattice) trellis.par.set(oc.theme()) d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d) tg <- tripGrid(tr) plot(tg)
oc.colors(10) library(lattice) trellis.par.set(oc.theme()) d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d) tg <- tripGrid(tr) plot(tg)
Trip rasterize.
x |
|
y |
Raster* object |
field |
attribute from which differences will be calculated, defaults to the time-stamp between trip locations |
RasterLayer
d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d, c("tms", "id")) tr$temp <- sort(runif(nrow(tr))) r <- rasterize(tr) rasterize(tr, grid = r) rasterize(tr, r, field = "temp") rasterize(tr, method = "density") rasterize(tr, method = "density", grid = r) rasterize(tr, r, field = "tms") rasterize(tr, r)
d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d, c("tms", "id")) tr$temp <- sort(runif(nrow(tr))) r <- rasterize(tr) rasterize(tr, grid = r) rasterize(tr, r, field = "temp") rasterize(tr, method = "density") rasterize(tr, method = "density", grid = r) rasterize(tr, r, field = "tms") rasterize(tr, r)
Return a (Spatial) data frame of location records from raw Argos files. Multiple files may be read, and each set of records is appended to the data frame in turn. Basic validation of the data is enforced by default.
readArgos( x, correct.all = TRUE, dtFormat = "%Y-%m-%d %H:%M:%S", tz = "GMT", duplicateTimes.eps = 0.01, p4 = "+proj=longlat +ellps=WGS84", verbose = FALSE, read_alt = NULL, ... ) readDiag(x, return_trip = FALSE, read_alt = 1L, ...)
readArgos( x, correct.all = TRUE, dtFormat = "%Y-%m-%d %H:%M:%S", tz = "GMT", duplicateTimes.eps = 0.01, p4 = "+proj=longlat +ellps=WGS84", verbose = FALSE, read_alt = NULL, ... ) readDiag(x, return_trip = FALSE, read_alt = 1L, ...)
x |
vector of file names of Argos "DAT" or "DIAG" files. |
correct.all |
logical - enforce validity of data as much as possible? (see Details) |
dtFormat |
the DateTime format used by the Argos data "date" and "time" pasted together |
tz |
timezone - GMT/UTC is assumed |
duplicateTimes.eps |
what is the tolerance for times being duplicate? |
p4 |
PROJ.4 projection string, "+proj=longlat +ellps=WGS84" is assumed |
verbose |
if TRUE, details on date-time adjustment is reported |
read_alt |
is |
... |
reserved for future use |
return_trip |
for |
readArgos
performs basic validation checks for class trip
are
made, and enforced based on correct.all
:
No duplicate records in the data, these are simply removed. Records are ordered by DateTime ("date", "time", "gmt") within ID ("ptt"). No duplicate DateTime values within ID are allowed: to enforce this the time values are moved forward by one second - this is done recursively and is not robust.
If validation fails the function will return a
SpatialPointsDataFrame-class
. Files that are not obviously of
the required format are skipped.
Argos location quality data "class" are ordered, assuming that the available
levels is levels=c("Z", "B", "A", "0", "1", "2", "3")
.
A projection string is added to the data, assuming the PROJ.4 longlat - if any longitudes are greater than 360 the PROJ.4 argument "+over" is added.
readDiag
simply builds a data.frame
.
With read_alt
the default value NULL
returns the PRV location as-is. Some files may have
a standardized location, and a dummy. If read_alt
is set to 1 or 2 the corresponding "alternative"
location is returned. 1 is a standardized location corresponding to the original PRV message, and 2 is a
"dummy" location.
readArgos
returns a trip
object, if all goes well, or simply a
SpatialPointsDataFrame-class
.
readDiag
returns a data.frame
with 8 columns:
lon1
,lat1
first pair of coordinates
lon1
,lat1
second pair of coordinates
gmt DateTimes as POSIXct
id Platform Transmitting Terminal (PTT) ID
lq Argos location quality class
iq some other thing
This works on some Argos files I have seen.
The Argos data documentation was (ca. 2003) at http://www.argos-system.org/manual. Specific details on the PRV ("provide data") format were found in Chapter 4_4_8, originally at 'http://www.cls.fr/manuel/html/chap4/chap4_4_8.htm'.
trip
, SpatialPointsDataFrame-class
,
adjust.duplicateTimes
, for manipulating these data, and
argos.sigma
for relating a numeric value to Argos quality
"classes".
sepIdGaps
for splitting the IDs in these data on some minimum
gap.
order
, duplicated
, , ordered
for
general manipulation of this type.
argosfile <- system.file("extdata/argos/98feb.dat", package = "trip", mustWork = TRUE) argos <- readArgos(argosfile)
argosfile <- system.file("extdata/argos/98feb.dat", package = "trip", mustWork = TRUE) argos <- readArgos(argosfile)
A reproj method for trip objects.
## S3 method for class 'trip' reproj(x, target, ..., source = NULL)
## S3 method for class 'trip' reproj(x, target, ..., source = NULL)
x |
trip object |
target |
target projection |
... |
ignored |
source |
projection of source data, usually ignore this for trips |
a trip reprojected to 'target'
Create a filter index of a track for "bad" points with a combination of speed, distance and angle tests.
sda(x, smax, ang = c(15, 25), distlim = c(2.5, 5), pre = NULL)
sda(x, smax, ang = c(15, 25), distlim = c(2.5, 5), pre = NULL)
x |
trip object |
smax |
maximum speed, in km/h |
ang |
minimum turning angle/s in degrees |
distlim |
maximum step lengths in km |
pre |
include this filter in the removal |
This is an independent implementation from that in the package argosfilter by Freitas 2008.
logical vector, with FALSE
values where the tests failed
Freitas, C., Lydersen, C., Fedak, M. A. and Kovacs, K. M. (2008), A simple new algorithm to filter marine mammal Argos locations. Marine Mammal Science, 24: 315?V325. doi: 10.1111/j.1748-7692.2007.00180.x
A new set of ID levels can be created by separating those given based on a minimum gap in another set of data. This is useful for separating instruments identified only by their ID into separate events in time.
sepIdGaps(id, gapdata, minGap = 3600 * 24 * 7)
sepIdGaps(id, gapdata, minGap = 3600 * 24 * 7)
id |
existing ID levels |
gapdata |
data matching |
minGap |
the minimum "gap" to use in gapdata to create a new ID level |
The assumption is that a week is a long time for a tag not to record anything.
A new set of ID levels, named following the pattern that "ID" split into 3 would provided "ID", "ID_2" and "ID_3".
It is assumed that each vector provides is sorted by gapdata
within
id
. No checking is done, and so it is suggested that this only be
used on ID columns within existing, validated trip
objects.
id <- gl(2, 8) gd <- Sys.time() + 1:16 gd[c(4:6, 12:16)] <- gd[c(4:6, 12:16)] + 10000 sepIdGaps(id, gd, 1000)
id <- gl(2, 8) gd <- Sys.time() + 1:16 gd[c(4:6, 12:16)] <- gd[c(4:6, 12:16)] + 10000 sepIdGaps(id, gd, 1000)
Create a filter of a track for "bad" points implying a speed of motion that is unrealistic.
speedfilter(x, max.speed = NULL, test = FALSE)
speedfilter(x, max.speed = NULL, test = FALSE)
x |
trip object |
max.speed |
speed in kilometres (or other unit) per hour, the unit is kilometres if the trip is in longitude latitude coordinates, or in the unit of the projection projection (usually metres per hour) |
test |
cut the algorithm short and just return first pass |
Using an algorithm (McConnnell et al., 1992), points are tested for speed between previous / next and 2nd previous / next points. Contiguous sections with an root mean square speed above a given maximum have their highest rms point removed, then rms is recalculated, until all points are below the maximum. By default an (internal) root mean square function is used, this can be specified by the user.
If the coordinates of the trip
data are not projected, or NA the
distance calculation assumes longlat and kilometres (great circle). For
projected coordinates the speed must match the units of the coordinate
system. (The PROJ.4 argument "units=km" is suggested).
Logical vector matching positions in the coordinate records that pass the filter.
This algorithm is destructive, and provides little information about location uncertainty. It is provided because it's commonly used and provides an illustrative benchmark for further work.
It is possible for the filter to become stuck in an infinite loop, depending on the function passed to the filter. Several minutes is probably too long for hundreds of points, test on smaller sections if unsure.
This algorithm was originally taken from IDL code by David Watts at the Australian Antarctic Division, and used in various other environments before the development of this version.
David Watts and Michael D. Sumner
The algorithm comes from McConnell, B. J. and Chambers, C. and Fedak, M. A. (1992) Foraging ecology of southern elephant seals in relation to the bathymetry and productivity of the southern ocean. Antarctic Science 4 393-398
sda
for a fast distance angle filter to combine with speed filtering
Object to identify DateTimes and IDs in a Spatial object.
TimeOrderedRecords(x)
TimeOrderedRecords(x)
x |
Character vector of 2 elements specifying the data columns of DateTimes and IDs |
TimeOrderedRecords
holds a 2-element character vector, naming the data columns
of DateTimes and IDs.
##' tor <- TimeOrderedRecords(c("datetime", "ID"))
##' tor <- TimeOrderedRecords(c("datetime", "ID"))
The main use of this class and creator function is for
SpatialPointsDataFrame-class
s which are used with
TimeOrderedRecords for the class trip
.
S4 object, TimeOrderedRecords (a class to hold the names of the date-time and id columns)
TOR.columns
:2-element vector of class "character"
Future versions may change significantly, this class is very basic and could probably be implemented in a better way. Specifying TOR columns by formula would be a useful addition.
TimeOrderedRecords
, trip
for creating trip objects, and trip-class
for that class
showClass("TimeOrderedRecords") tor <- new("TimeOrderedRecords", TOR.columns=c("datetime", "ID"))
showClass("TimeOrderedRecords") tor <- new("TimeOrderedRecords", TOR.columns=c("datetime", "ID"))
Calculate the angles between subsequent 2-D coordinates using Great Circle distance (spherical) methods.
trackAngle(x) ## S3 method for class 'trip' trackAngle(x) ## Default S3 method: trackAngle(x)
trackAngle(x) ## S3 method for class 'trip' trackAngle(x) ## Default S3 method: trackAngle(x)
x |
trip object, or matrix of 2-columns, with x/y coordinates |
If x
is a trip object, the return result has an extra element for the
start and end point of each individual trip, with value NA.
This is an optimized hybrid of "raster::bearing" and "maptools::gzAzimuth". New code is in the traipse package.
Vector of angles (degrees) between coordinates.
Calculate the distances between subsequent 2-D coordinates using Euclidean or Great Circle distance (WGS84 ellipsoid) methods.
trackDistance(x1, y1, x2, y2, longlat = TRUE, prev = FALSE)
trackDistance(x1, y1, x2, y2, longlat = TRUE, prev = FALSE)
x1 |
trip object, matrix of 2-columns, with x/y coordinates OR a vector of x start coordinates |
y1 |
vector of y start coordinates, if x1 is not a matrix |
x2 |
vector of x end coordinates, if x1 is not a matrix |
y2 |
vector of y end coordinates, if x1 is not a matrix |
longlat |
if FALSE, Euclidean distance, if TRUE Great Circle distance |
prev |
if TRUE and x1 is a trip, the return value has a padded end value (\"prev\"ious), rather than start (\"next\") |
If x1
is a trip object, arguments x2
, x3
, y2
are
ignored and the return result has an extra element for the start point of
each individual trip, with value 0.0.
The prev
argument is ignore unless x1 is a trip.
Distance values are in the units of the input coordinate system when longlat is FALSE, and in kilometres when longlat is TRUE.
This originally used spDistsN1
, then implemented the sp
gcdist
source directly in R, and now uses geodist
.
Please see the traipse package for a more modern approach.
Vector of distances between coordinates.
Original source taken from sp package, but now using Helmert from Karney (2013) see the geodist package.
d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d, c("tms", "id")) ## the method knows this is a trip, so there is a distance for every ## point, including 0s as the start and at transitions between ## individual trips trackDistance(tr) ## the default method does not know about the trips, so this is ##(n-1) distances between all points trackDistance(coordinates(tr), longlat = FALSE) ## we get NA at the start, end and at transitions between trips angles <- trackAngle(tr)
d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d, c("tms", "id")) ## the method knows this is a trip, so there is a distance for every ## point, including 0s as the start and at transitions between ## individual trips trackDistance(tr) ## the default method does not know about the trips, so this is ##(n-1) distances between all points trackDistance(coordinates(tr), longlat = FALSE) ## we get NA at the start, end and at transitions between trips angles <- trackAngle(tr)
Functions for retrieving the names of the columns used for DateTime and ID, as well as the data.
getTORnames(obj) getTimeID(obj) ## S3 method for class 'summary.TORdata' print(x, ...)
getTORnames(obj) getTimeID(obj) ## S3 method for class 'summary.TORdata' print(x, ...)
obj |
|
x |
trip object |
... |
currently ignored |
getTORnames
retrieves the column names from an object extending the
class TimeOrderedRecords
, and getTimeID
returns the data as a
data frame from an object extending the class TimeOrderedRecords
.
trip-class
, for the use of this class with
SpatialPointsDataFrame-class
.
tor <- TimeOrderedRecords(c("time", "id")) getTORnames(tor)
tor <- TimeOrderedRecords(c("time", "id")) getTORnames(tor)
An extension of SpatialPointsDataFrame-class
by including
"TimeOrderedRecords"
. The records within the data frame are
explicitly ordered by DateTime data within IDs.
Objects can be created by calls of the form
trip(obj="SpatialPointsDataFrame", TORnames="TimeOrderedRecords")
.
The object contains all the slots present within a
SpatialPointsDataFrame-class
, particularly data
which
contains columns of at least those specified by TOR.columns
.
trip
for examples of directly using the class.
trip-accessors
describes methods for accessing information on
trip
objects.
showClass("trip") d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d) summary(tr) plot(tr) lines(tr) dim(tr) names(tr) subset(tr, id == "2") as.data.frame(tr) tr[1:3, ] tr[, 1] tr[[1]]
showClass("trip") d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) tr <- trip(d) summary(tr) plot(tr) lines(tr) dim(tr) names(tr) subset(tr, id == "2") as.data.frame(tr) tr[1:3, ] tr[, 1] tr[[1]]
trip
objectsCreate an object of class trip
, extending the basic functionality
of SpatialPointsDataFrame-class
by specifying the data columns
that define the "TimeOrdered" quality of the records.
trip(obj, TORnames, correct_all = TRUE) trip(obj) <- value ## S4 method for signature 'trip,ANY' split(x, f, drop = FALSE, ...) ## S4 method for signature 'trip,ANY,ANY,ANY' x[i, j, ..., drop = TRUE]
trip(obj, TORnames, correct_all = TRUE) trip(obj) <- value ## S4 method for signature 'trip,ANY' split(x, f, drop = FALSE, ...) ## S4 method for signature 'trip,ANY,ANY,ANY' x[i, j, ..., drop = TRUE]
obj |
A data frame, a grouped data frame or a |
TORnames |
Either a |
correct_all |
logical value, if |
value |
A 4-element character vector specifying the X, Y, DateTime coordinates
and ID of |
x |
trip object |
f |
grouping vector as per |
drop |
unused but necessary for method consistency |
i , j , ...
|
indices specifying elements to extract |
The original form of trip()
required very strict input as a 'SpatialPointsDataFrame' and
specifying which were the time and ID columns, but the input can be more flexible. If the object is a
grouped data frame ('dplyr-style') then the (first) grouping is assumed to define individual trips and that
columns 1, 2, 3 are the x-, y-, time-coordinates in that order. It can also be a trip
object for
redefining TORnames
.
The trip()
function can ingest track_xyt
, telemetry
, SpatialPointsDataFrame
, sf
,
trackeRdata
, grouped_df
, data.frame
, tbl_df
, mousetrap
, and in some cases
lists of those objects. Please get in touch if you think something that should work does not.
Track data often contains problems, with missing values in location or time,
times out of order or with duplicated times. The correct_all
argument is
set to TRUE
by default and will report any inconsistencies. Data really should
be checked first rather than relying on this auto-cleanup. The following problems are common:
duplicated records (every column with the same value in another row)
duplicated date-time values
missing date-time values, or missing x or y coordinates
records out of order within trip ID
For some data types there's no formal structure, but a simple convention such as
a set of names in a data frame. For example, the VTrack package has AATAMS1
which may be
turned into a trip with
trip(AATAMS1 %>% dplyr::select(longitude, latitude, timestamp, tag.ID, everything())
In time we can add support for all kinds of variants, detected by the names and contents.
See Chapter 2 of the trip thesis for more details.
A trip object, with the usual slots of a
SpatialPointsDataFrame-class
and the added
TimeOrderedRecords
. For the most part this can be treated as a
data.frame
with Spatial
coordinates.
Most of the methods available are by virtue of the sp package. Some, such
as split.data.frame
have been added to SPDF so that trip has the same
functionality.
signature(obj="SpatialPointsDataFrame",
TORnames="ANY")
The main construction.
signature(obj="SpatialPointsDataFrame",
TORnames="TimeOrderedRecords")
Object and TimeOrdered records class
signature(obj="ANY", TORnames="TimeOrderedRecords")
:
create a trip
object from a data frame.
signature(obj="trip", TORnames="ANY")
: (Re)-create a
trip
object using a character vector for TORnames
.
signature(obj="trip", TORnames="TimeOrderedRecords")
:
(re)-create a trip object using a TimeOrderedRecords
object.
speedfilter
, and tripGrid
for simplistic
speed filtering and spatial time spent gridding.
d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) ## the simplest way to create a trip is by order of columns trip(d) tr <- trip(d) ## real world data in CSV mi_dat <- read.csv(system.file("extdata/MI_albatross_sub10.csv", package = "trip"), stringsAsFactors = FALSE) mi_dat$gmt <- as.POSIXct(mi_dat$gmt, tz = "UTC") mi_dat$sp_id <- sprintf("%s%s_%s_%s", mi_dat$species, substr(mi_dat$breeding_status, 1, 1), mi_dat$band, mi_dat$tag_ID) sp::coordinates(mi_dat) <- c("lon", "lat") ## there are many warnings, but the outcome is fine ## (sp_id == 'WAi_14030938_2123' has < 3 locations as does LMi_12143650_14257) mi_dat <- trip(mi_dat, c("gmt", "sp_id") ) plot(mi_dat, pch = ".") #lines(mi_dat) ## ugly mi_dat_polar <- reproj(mi_dat, "+proj=stere +lat_0=-90 +lon_0=154 +datum=WGS84") plot(mi_dat_polar, pch = ".") lines(mi_dat_polar)
d <- data.frame(x=1:10, y=rnorm(10), tms=Sys.time() + 1:10, id=gl(2, 5)) ## the simplest way to create a trip is by order of columns trip(d) tr <- trip(d) ## real world data in CSV mi_dat <- read.csv(system.file("extdata/MI_albatross_sub10.csv", package = "trip"), stringsAsFactors = FALSE) mi_dat$gmt <- as.POSIXct(mi_dat$gmt, tz = "UTC") mi_dat$sp_id <- sprintf("%s%s_%s_%s", mi_dat$species, substr(mi_dat$breeding_status, 1, 1), mi_dat$band, mi_dat$tag_ID) sp::coordinates(mi_dat) <- c("lon", "lat") ## there are many warnings, but the outcome is fine ## (sp_id == 'WAi_14030938_2123' has < 3 locations as does LMi_12143650_14257) mi_dat <- trip(mi_dat, c("gmt", "sp_id") ) plot(mi_dat, pch = ".") #lines(mi_dat) ## ugly mi_dat_polar <- reproj(mi_dat, "+proj=stere +lat_0=-90 +lon_0=154 +datum=WGS84") plot(mi_dat_polar, pch = ".") lines(mi_dat_polar)
These functions will be declared defunct in a future release.
as.SpatialLinesDataFrame.trip(from) trip.split.exact(x, dates) as.ltraj.trip(xy) as.trip.SpatialLinesDataFrame(from)
as.SpatialLinesDataFrame.trip(from) trip.split.exact(x, dates) as.ltraj.trip(xy) as.trip.SpatialLinesDataFrame(from)
from |
trip object |
x |
see |
dates |
see |
xy |
|
Create a grid of time spent from an object of class trip
by exact
cell crossing methods, weighted by the time between locations for separate
trip events.
tripGrid(x, grid = NULL, method = "pixellate", ...)
tripGrid(x, grid = NULL, method = "pixellate", ...)
x |
object of class |
grid |
GridTopology - will be generated automatically if NULL |
method |
pixellate or density |
... |
pass arguments to density.psp if that method is chosen (and
temporary mechanism to direct users of legacy methods to
|
Zero-length lines cannot be summed directly, their time value is summed by assuming the line is a point. A warning used to be given, but as it achieved nothing but create confusion it has been removed. The density method returns proportionate values, not summed time durations.
See pixellate.psp
and pixellate.ppp
for the details on the
method used. See density.psp
for method="density".
Trip events are assumed to start and end as per the object passed in. To
work with inferred "cutoff" positions see split.trip.exact
.
tripGrid
returns an object of class SpatialGridDataFrame
, with
one column "z" containing the time spent in each cell in seconds.
Create a grid of time spent from an object of class trip
by
approximating the time between locations for separate trip events.
tripGrid.interp(x, grid = NULL, method = "count", dur = NULL, ...) kdePoints(x, h = NULL, grid = NULL, resetTime = TRUE, ...) countPoints(x, dur = 1, grid = NULL)
tripGrid.interp(x, grid = NULL, method = "count", dur = NULL, ...) kdePoints(x, h = NULL, grid = NULL, resetTime = TRUE, ...) countPoints(x, dur = 1, grid = NULL)
x |
object of class trip |
grid |
GridTopology - will be generated automatically if NULL |
method |
name of method for quantifying time spent, see Details |
dur |
The \"dur\"ation of time used to interpolate between available locations (see Details) |
... |
other arguments passed to |
h |
kernel bandwidth |
resetTime |
rescale result back to the total duration of the input |
This set of functions was the the original tripGrid from prior to version
1.1-6. tripGrid
should be used for more exact and fast calculations
assuming linear motion between fixes.
The intention is for tripGrid.interp
to be used for exploring
approximate methods of line-to-cell gridding.
Trip locations are first interpolated, based on an equal-time spacing
between records. These interpolated points are then "binned" to a grid of
cells. The time spacing is specified by the dur
(duration) argument to
interpequal
in seconds (i.e. dur=3600
is used for 1 hour).
Shorter time periods will require longer computation with a closer
approximation to the total time spent in the gridded result.
Currently there are methods "count" and "kde" for quantifying time spent, corresponding to the functions "countPoints" and "kdePoints". "kde" uses kernel density to smooth the locations, "count" simply counts the points falling in a grid cell.
tripGrid
returns an object of class SpatialGridDataFrame
, with
one column "z" containing the time spent in each cell in seconds. If
kdePoints is used the units are not related to the time values and must be
scaled for further use.
bandwidth.nrd
for the calculation of bandwidth values used internally when not supplied by the user
Behavior of Pacific Walruses Tracked from the Alaska Coast of the Chukchi Sea.
Data set is provided as a 'trip' object. This is the abstract for the work:
"We tracked movements and haulout foraging behavior of walruses instrumented with satellite-linked data loggers from the Alaskan shores of the Chukchi Sea during the autumn of 2009 (n=13) and 2010 (n=2)." Jay, C. V. and Fischbach, A.S.
data(walrus818) plot(walrus818) lines(walrus818)
data(walrus818) plot(walrus818) lines(walrus818)
A spatial polygons object with coastlines of the northern hemisphere.
world_north
world_north
An object of class SpatialPolygonsDataFrame
with 185 rows and 11 columns.
This data set exists purely to avoid requiring reprojection in the vignette, the data uses the same projection as walrus818.
Export track data to a KML file, for use in Google Earth the continuous time slider.
write_track_kml( id, lon, lat, utc, z = NULL, kml_file = tempfile(fileext = ".kmz"), name = NULL, altitude_mode = c("absolute", "clampToGround", "clampToSeaFloor", "relativeToGround", "relativeToSeaFloor") )
write_track_kml( id, lon, lat, utc, z = NULL, kml_file = tempfile(fileext = ".kmz"), name = NULL, altitude_mode = c("absolute", "clampToGround", "clampToSeaFloor", "relativeToGround", "relativeToSeaFloor") )
id |
vector of grouping IDs (or a trip object) |
lon |
vector of longitude (ignored if id is a trip) |
lat |
vector of latitude (ignored if id is a trip) |
utc |
vector of POSIXct date-times (ignored if id is a trip) |
z |
vector of elevations, this cannot be set if 'id' is a trip |
kml_file |
filename for KML (KML or KMZ) (must end in .kml or .kmz) |
name |
internal name of dat (derived from kml_file if not specified) |
altitude_mode |
the altitude mode, 'absolute', 'clampToGround', 'clampToSeaFloor', 'relativeToGround', or 'relativeToSeaFloor', see Details |
To include altitude set every argument explicitly, by input of separate 'id', 'lon', 'lat', 'utc' and 'z' arguments. If the first argument 'id' is a trip object there is no facility to include the 'z' altitude values.
If 'z' is included it is applied as a third coordinate, with 'altitude_mode' controlling the interpretation, see https://developers.google.com/kml/documentation/altitudemode. If the 'kml_file' ends with ".kmz" the file is compressed, otherwise it must end with ".kml" and the compression archive step is not applied.
Sadly the interactive time slider is only available with the desktop version of Google Earth, the data loads into the browser version but can't be interactive.
character vector, file name location of file produced
Original implementation by Tomislav Hengl in the 'plotKML' package for 'SpatialLinesDataFrame', adapted by M. Sumner for use in continuous-time form.
kfile <- write_track_kml(walrus818[seq(1, 1000, by = 5), ]) print(kfile) unlink(kfile)
kfile <- write_track_kml(walrus818[seq(1, 1000, by = 5), ]) print(kfile) unlink(kfile)