Title: | Empirically Informed Random Trajectory Generation in 3-D |
---|---|
Description: | Creates realistic random trajectories in a 3-D space between two given fix points, so-called conditional empirical random walks (CERWs). The trajectory generation is based on empirical distribution functions extracted from observed trajectories (training data) and thus reflects the geometrical movement characteristics of the mover. A digital elevation model (DEM), representing the Earth's surface, and a background layer of probabilities (e.g. food sources, uplift potential, waterbodies, etc.) can be used to influence the trajectories. Unterfinger M (2018). "3-D Trajectory Simulation in Movement Ecology: Conditional Empirical Random Walk". Master's thesis, University of Zurich. <https://www.geo.uzh.ch/dam/jcr:6194e41e-055c-4635-9807-53c5a54a3be7/MasterThesis_Unterfinger_2018.pdf>. Technitis G, Weibel R, Kranstauber B, Safi K (2016). "An algorithm for empirically informed random trajectory generation between two endpoints". GIScience 2016: Ninth International Conference on Geographic Information Science, 9, online. <doi:10.5167/uzh-130652>. |
Authors: | Merlin Unterfinger [aut, cre] |
Maintainer: | Merlin Unterfinger <[email protected]> |
License: | GPL-3 |
Version: | 0.7.0 |
Built: | 2025-02-01 04:14:09 UTC |
Source: | https://github.com/munterfi/ertg3d |
Calculates the chi maps for one rasterStack
or all raster all the raster pairs stored in two rasterStack
s.
As observed values, the first stack is used. The expected value is either set to the mean of the first stack, or if given
to be the values of the second stack.
chiMaps(stack1, stack2 = NULL, verbose = FALSE)
chiMaps(stack1, stack2 = NULL, verbose = FALSE)
stack1 |
|
stack2 |
|
verbose |
logical: print currently processed height band in raster stack? |
A rasterStack
containing the chi maps.
print("tbd.")
print("tbd.")
This is data to be included in the package and can be used to test its functionality.
The dem
data is a RasterLayer
and has a resolution of 90 meters. It is the topography
of the Swiss midlands. The complete dataset can be downloaded directly from http://srtm.csi.cgiar.org/srtmdata/.
https://srtm.csi.cgiar.org/srtmdata/
Crops the DEM to the extent of the track with a buffer
dem2track.extent(DEM, track, buffer = 100)
dem2track.extent(DEM, track, buffer = 100)
DEM |
a raster containing a digital elevation model, covering the extent as the track |
track |
data.frame with x,y,z coordinates of the original track |
buffer |
buffer with, by default set to 100 |
A the cropped digital elevation model as a raster layer.
dem2track.extent(dem, niclas)
dem2track.extent(dem, niclas)
Distance of each track point to a given point
dist2point.3d(track, point, groundDistance = FALSE)
dist2point.3d(track, point, groundDistance = FALSE)
track |
a list containing data.frames with x,y,z coordinates or a data.frame |
point |
a vector with x, y or x, y, z coordinates |
groundDistance |
logical: calculate only ground distance in x-y plane? |
Returns the distance of each track point to the point.
dist2point.3d(niclas, c(0, 0, 0))
dist2point.3d(niclas, c(0, 0, 0))
Calculates the distance between every point in the track and the last point (target).
dist2target.3d(track)
dist2target.3d(track)
track |
a track data.frame containing x, y and z coordinates |
A numeric vector with the distances to target
dist2target.3d(niclas)
dist2target.3d(niclas)
The empirically informed random trajectory generator in three dimensions (eRTG3D) is an algorithm to generate realistic random trajectories in a 3-D space between two given fix points in space, so-called Conditional Empirical Random Walks. The trajectory generation is based on empirical distribution functions extracted from observed trajectories (training data) and thus reflects the geometrical movement characteristics of the mover. A digital elevation model (DEM), representing the Earth's surface, and a background layer of probabilities (e.g. food sources, uplift potential, waterbodies, etc.) can be used to influence the trajectories.
See the packages site on GitHub, detailed information about the algorithm in this Master’s Thesis, or test the algorithm online in the eRTG3D Simulator.
Function to filter out tracks that have found a dead end
filter.dead.ends(cerwList)
filter.dead.ends(cerwList)
cerwList |
list of data.frames and NULL entries |
A list that is only containing valid tracks.
filter.dead.ends(list(niclas, niclas))
filter.dead.ends(list(niclas, niclas))
Creates a list consisting of the three dimensional probability distribution cube for turning angle, lift angle and step length (turnLiftStepHist) as well as the uni-dimensional distributions of the differences of the turn angles, lift angles and step lengths with a lag of 1 to maintain minimal level of autocorrelation in each of the terms. Additionally also the distribution of the flight height over the ellipsoid (absolute) and the distribution of flight height over the topography (relative) can be included.
get.densities.3d( turnAngle, liftAngle, stepLength, deltaLift, deltaTurn, deltaStep, gradientAngle = NULL, heightEllipsoid = NULL, heightTopo = NULL, maxBin = 25 )
get.densities.3d( turnAngle, liftAngle, stepLength, deltaLift, deltaTurn, deltaStep, gradientAngle = NULL, heightEllipsoid = NULL, heightTopo = NULL, maxBin = 25 )
turnAngle |
turn angles of the track (t) |
liftAngle |
lift angles of the track (l) |
stepLength |
stepLength of the track (d) |
deltaLift |
auto differences of the turn angles (diff(t)) |
deltaTurn |
auto differences of the lift angles (diff(l)) |
deltaStep |
auto differences of the step length (diff(d)) |
gradientAngle |
|
heightEllipsoid |
flight height over the ellipsoid (absolute) or |
heightTopo |
flight height over the topography (relative) or |
maxBin |
numeric scalar, maximum number of bins per dimension of the tld-cube (turnLiftStepHist) |
A list containing the tldCube and the autodifferences functions (and additionally the flight height distribution functions)
niclas <- track.properties.3d(niclas)[2:nrow(niclas), ] P <- get.densities.3d( turnAngle = niclas$t, liftAngle = niclas$l, stepLength = niclas$d, deltaLift = diff(niclas$t), deltaTurn = diff(niclas$l), deltaStep = diff(niclas$d), gradientAngle = NULL, heightEllipsoid = NULL, heightTopo = NULL, maxBin = 25 )
niclas <- track.properties.3d(niclas)[2:nrow(niclas), ] P <- get.densities.3d( turnAngle = niclas$t, liftAngle = niclas$l, stepLength = niclas$d, deltaLift = diff(niclas$t), deltaTurn = diff(niclas$l), deltaStep = diff(niclas$d), gradientAngle = NULL, heightEllipsoid = NULL, heightTopo = NULL, maxBin = 25 )
Calculates the ratio between horizontal movement and vertical movement. The value expresses the distance covered forward movement per distance movement in sinking.
get.glideRatio.3d(track)
get.glideRatio.3d(track)
track |
a track data.frame containing x, y and z coordinates of a gliding section |
The ratio between horizontal and vertical movement.
get.glideRatio.3d(niclas)
get.glideRatio.3d(niclas)
Creates a list consisting of the 3 dimensional probability distribution cube for turning angle, lift angle and step length (turnLiftStepHist) as well as the uni-dimensional distributions of the differences of the turning angles, lift angles and step lengths with a lag of 1 to maintain minimal level of autocorrelation in each of the terms.
get.section.densities.3d( trackSections, gradientDensity = TRUE, heightDistEllipsoid = TRUE, DEM = NULL, maxBin = 25 )
get.section.densities.3d( trackSections, gradientDensity = TRUE, heightDistEllipsoid = TRUE, DEM = NULL, maxBin = 25 )
trackSections |
list of track sections got by the track.split.3d function |
gradientDensity |
logical: Should a distribution of the gradient angle be extracted and later used in the simulations? |
heightDistEllipsoid |
logical: Should a distribution of the flight height over ellipsoid be extracted and later used in the sim.cond.3d()? |
DEM |
a raster containing a digital elevation model, covering the same extent as the track sections |
maxBin |
numeric scalar, maximum number of bins per dimension of the tld-cube (turnLiftStepHist) |
A list containing the tldCube and the autodifferences functions (and additionally the height distribution function)
get.section.densities.3d(list(niclas[1:10, ], niclas[11:nrow(niclas), ]))
get.section.densities.3d(list(niclas[1:10, ], niclas[11:nrow(niclas), ]))
Get densities creates a list consisting of the 3 dimensional probability distribution cube for turning angle, lift angle and step length (turnLiftStepHist) as well as the uni-dimensional distributions of the differences of the turning angles, lift angles and step lengths with a lag of 1 to maintain minimal level of autocorrelation in each of the terms.
get.track.densities.3d( track, gradientDensity = TRUE, heightDistEllipsoid = TRUE, DEM = NULL, maxBin = 25 )
get.track.densities.3d( track, gradientDensity = TRUE, heightDistEllipsoid = TRUE, DEM = NULL, maxBin = 25 )
track |
a data.frame with 3 columns containing the x,y,z coordinates |
gradientDensity |
logical: Should a distribution of the gradient angle be extracted and later used in the simulations? |
heightDistEllipsoid |
logical: Should a distribution of the flight height over ellipsoid be extracted and later used in the sim.cond.3d()? |
DEM |
a raster containing a digital elevation model, covering the same extent as the track |
maxBin |
numeric scalar, maximum number of bins per dimension of the tld-cube (turnLiftStepHist) |
A list containing the tldCube and the autodifferences functions (and additionally the height distribution function)
The time between the acquisition of fix points of the track must be constant, otherwise this leads to distorted statistic distributions, which increases the probability of dead ends. In this case please check track.split.3d and get.section.densities.3d
get.track.densities.3d(niclas, heightDist = TRUE)
get.track.densities.3d(niclas, heightDist = TRUE)
'sf, data.frame'
)Tests if the object is a simple feature collection (class: 'sf, data.frame'
)
is.sf.3d(track)
is.sf.3d(track)
track |
any object to test |
A logical: TRUE
if is a simple feature collection (class: 'sf, data.frame'
) of the sf package, FALSE
otherwise.
is.sf.3d(niclas) is.sf.3d(track2sf.3d(track = niclas, CRS = 2056))
is.sf.3d(niclas) is.sf.3d(track2sf.3d(track = niclas, CRS = 2056))
Calculates the lift angle between every point in the track and the last point (target).
lift2target.3d(track)
lift2target.3d(track)
track |
a track data.frame containing x, y and z coordinates |
A numeric vector with the lift angles to target
lift2target.3d(niclas)
lift2target.3d(niclas)
Avoids the problem of -Inf occuring for log(0).
logRasterStack(rStack, standartize = FALSE, InfVal = NA)
logRasterStack(rStack, standartize = FALSE, InfVal = NA)
rStack |
rasterStack to convert to logarithmic scale |
standartize |
logical: standartize cube between 0 and 1 |
InfVal |
the value that |
A rasterStack in logarithmic scale
logRasterStack(raster::stack(dem))
logRasterStack(raster::stack(dem))
Applies a twosided moving median window on a vector, where the window paramter is the total size of the window. The value in the window middle is the index where the median of the window is written. Therefore the window size has to be an uneven number. The border region of the vetor is filled with a one-sided median. There might be border effects.
movingMedian(data, window)
movingMedian(data, window)
data |
numeric vector |
window |
uneven number for the size of the moving window |
A numeric vector.
movingMedian(sequence(1:10), window = 5)
movingMedian(sequence(1:10), window = 5)
Creates multiple conditional empirical random walks, with a specific starting and ending point, geometrically similar to the initial trajectory by applying sim.cond.3d multiple times.
n.sim.cond.3d( n.sim, n.locs, start = c(0, 0, 0), end = start, a0, g0, densities, qProbs, error = FALSE, parallel = FALSE, DEM = NULL, BG = NULL )
n.sim.cond.3d( n.sim, n.locs, start = c(0, 0, 0), end = start, a0, g0, densities, qProbs, error = FALSE, parallel = FALSE, DEM = NULL, BG = NULL )
n.sim |
number of CERWs to simulate |
n.locs |
length of the trajectory in locations |
start |
numeric vector of length 3 with the coordinates of the start point |
end |
numeric vector of length 3 with the coordinates of the end point |
a0 |
initial incoming heading in radian |
g0 |
initial incoming gradient/polar angle in radian |
densities |
list object returned by the get.densities.3d function |
qProbs |
list object returned by the qProb.3d function |
error |
logical: add random noise to the turn angle, lift angle and step length to account for errors measurements? |
parallel |
logical: run computations in parallel (n-1 cores)? Or numeric: the number of nodes (maximum: n - 1 cores) |
DEM |
raster layer containing a digital elevation model, covering the area between start and end point |
BG |
a background raster layer that can be used to inform the choice of steps |
A list containing the CERWs or NULL
s if dead ends have been encountered.
niclas <- track.properties.3d(niclas) n.locs <- 3 P <- get.track.densities.3d(niclas) f <- 1500 start <- Reduce(c, niclas[1, 1:3]) end <- Reduce(c, niclas[n.locs, 1:3]) a0 <- niclas$a[1] g0 <- niclas$g[1] uerw <- sim.uncond.3d( n.locs * f, start = start, a0 = a0, g0 = g0, densities = P ) Q <- qProb.3d(uerw, n.locs) n.sim.cond.3d( n.sim = 2, n.locs = n.locs, start = start, end = end, a0 = a0, g0 = g0, densities = P, qProbs = Q )
niclas <- track.properties.3d(niclas) n.locs <- 3 P <- get.track.densities.3d(niclas) f <- 1500 start <- Reduce(c, niclas[1, 1:3]) end <- Reduce(c, niclas[n.locs, 1:3]) a0 <- niclas$a[1] g0 <- niclas$g[1] uerw <- sim.uncond.3d( n.locs * f, start = start, a0 = a0, g0 = g0, densities = P ) Q <- qProb.3d(uerw, n.locs) n.sim.cond.3d( n.sim = 2, n.locs = n.locs, start = start, end = end, a0 = a0, g0 = g0, densities = P, qProbs = Q )
Creates conditional empirical random walks in gliding mode, between a start and end point. The walk is performed on a MODE layer and, if provided, additionally on a background and digital elevation layer. The gliding is simulated with sim.cond.3d and soaring with sim.uncond.3d, therefore soaring is not restricted towards the target and can happen completly free as long as there are good thermal conditions. It is important to extract for every mode in the MODE raster layer a corresponding densities object with get.densities.3d and pass them to the function.
n.sim.glidingSoaring.3d( n.sim = 1, parallel = FALSE, MODE, dGliding, dSoaring, qGliding, start = c(0, 0, 0), end = start, a0, g0, error = TRUE, smoothTransition = TRUE, glideRatio = 20, DEM = NULL, BG = NULL, verbose = FALSE )
n.sim.glidingSoaring.3d( n.sim = 1, parallel = FALSE, MODE, dGliding, dSoaring, qGliding, start = c(0, 0, 0), end = start, a0, g0, error = TRUE, smoothTransition = TRUE, glideRatio = 20, DEM = NULL, BG = NULL, verbose = FALSE )
n.sim |
number of simulations to produce |
parallel |
logical: run computations in parallel (n-1 cores)? Or numeric: the number of nodes (maximum: n - 1 cores) |
MODE |
raster layer containing the number/index of the mode, which should be used at each location |
dGliding |
density object returned by the get.densities.3d function for gliding mode |
dSoaring |
density object returned by the get.densities.3d function for soaring mode |
qGliding |
the Q probabilites for the steps in gliding mode (qProb.3d) |
start |
numeric vector of length 3 with the coordinates of the start point |
end |
numeric vector of length 3 with the coordinates of the end point |
a0 |
initial incoming heading in radian |
g0 |
initial incoming gradient/polar angle in radian |
error |
logical: add random noise to the turn angle, lift angle and step length to account for errors measurements? |
smoothTransition |
logical: should the transitions between soaring and the following gliding sections be smoothed? Recommended to avoid dead ends |
glideRatio |
ratio between vertical and horizontal movement, by default set to 15 meters forward movement per meter vertical movement |
DEM |
raster layer containing a digital elevation model, covering the area between start and end point |
BG |
a background raster layer that can be used to inform the choice of steps |
verbose |
logical: print current mode used? |
A list containing 'soaring-gliding' trajectories or NULL
s if dead ends have been encountered.
The MODE raster layer must be in the following structure: Gliding pixels have the value 1 and soaring pixel the values 2. NA
's are not allowed in the raster.
print("tbd.")
print("tbd.")
This is data to be included in the package and can be used to test its functionality.
The track consists of x, y and z coordinates and represents the movement of a stork
called niclas
in the Swiss midlands.
Function detects the operating system and chooses the approximate kind of process for parallelizing the task: Windows: PSOCKCluster, Unix: Forking.
parpbapply( X, FUN, MARGIN, packages = NULL, export = NULL, envir = environment(), nNodes = parallel::detectCores() - 1 )
parpbapply( X, FUN, MARGIN, packages = NULL, export = NULL, envir = environment(), nNodes = parallel::detectCores() - 1 )
X |
an array, including a matrix. |
FUN |
function, the function to be applied to each element of X |
MARGIN |
a vector giving the subscripts which the function will be applied over. E.g., for a matrix 1 indicates rows, 2 indicates columns, c(1, 2) indicates rows and columns. Where X has named dimnames, it can be a character vector selecting dimension names. |
packages |
character vector, Only relevant for Windows: the packages needed in the function provided, eg. c("MASS", "data.table") |
export |
character vector, Only relevant for Windows: the varibales needed in the function provided, eg. c("df", "vec") |
envir |
environment, Only relevant for Windows: Environment from which the variables should be exported from |
nNodes |
numeric, Number of processes to start (unix: best to fit with the available Cores) |
Returns a vector or array or list of values obtained by applying a function to margins of an array or matrix.
n <- 1000 df <- data.frame( x = seq(1, n, 1), y = -seq(1, n, 1) ) a <- parpbapply(X = df, FUN = sum, MARGIN = 1, nNodes = 2)
n <- 1000 df <- data.frame( x = seq(1, n, 1), y = -seq(1, n, 1) ) a <- parpbapply(X = df, FUN = sum, MARGIN = 1, nNodes = 2)
Function detects the operating system and chooses the approximate kind of process for parallelizing the task: Windows: PSOCKCluster, Unix: Forking.
parpblapply( X, FUN, packages = NULL, export = NULL, envir = environment(), nNodes = parallel::detectCores() - 1 )
parpblapply( X, FUN, packages = NULL, export = NULL, envir = environment(), nNodes = parallel::detectCores() - 1 )
X |
a vector (atomic or list) or an expression object. Other objects (including classed objects) will be coerced by base::as.list |
FUN |
function, the function to be applied to each element of X |
packages |
character vector, Only relevant for Windows: the packages needed in the function provided, eg. c("MASS", "data.table") |
export |
character vector, Only relevant for Windows: the varibales needed in the function provided, eg. c("df", "vec") |
envir |
environment, Only relevant for Windows: Environment from which the variables should be exported from |
nNodes |
numeric, Number of processes to start (unix: best to fit with the available Cores) |
A list with the results.
square <- function(x) { x * x } l <- parpblapply(X = 1:1000, FUN = square, export = c("square"), nNodes = 2)
square <- function(x) { x * x } l <- parpblapply(X = 1:1000, FUN = square, export = c("square"), nNodes = 2)
Function detects the operating system and chooses the approximate kind of process for parallelizing the task: Windows: PSOCKCluster, Unix: Forking.
parpbsapply( X, FUN, packages = NULL, export = NULL, envir = environment(), nNodes = parallel::detectCores() - 1 )
parpbsapply( X, FUN, packages = NULL, export = NULL, envir = environment(), nNodes = parallel::detectCores() - 1 )
X |
a vector (atomic or list) or an expression object. Other objects (including classed objects) will be coerced by base::as.list. |
FUN |
function, the function to be applied to each element of X |
packages |
character vector, Only relevant for Windows: the packages needed in the function provided, eg. c("MASS", "data.table") |
export |
character vector, Only relevant for Windows: the varibales needed in the function provided, eg. c("df", "vec") |
envir |
environment, Only relevant for Windows: Environment from which the variables should be exported from |
nNodes |
numeric, Number of processes to start (unix: best to fit with the available Cores) |
A vector with the results.
square <- function(x) { x * x } s <- parpbsapply(X = 1:1000, FUN = square, export = c("square"), nNodes = 2)
square <- function(x) { x * x } s <- parpbsapply(X = 1:1000, FUN = square, export = c("square"), nNodes = 2)
Plot function to plot the 3-D tracks in 2-D plane
plot2d( origTrack, simTrack = NULL, titleText = character(1), DEM = NULL, BG = NULL, padding = 0.1, alpha = 0.7, resolution = 500 )
plot2d( origTrack, simTrack = NULL, titleText = character(1), DEM = NULL, BG = NULL, padding = 0.1, alpha = 0.7, resolution = 500 )
origTrack |
a list containing data.frames with x,y,z coordinates or a data.frame |
simTrack |
a list containing data.frames with x,y,z coordinates or a data.frame |
titleText |
string with title of the plot |
DEM |
an object of type |
BG |
an object of type |
padding |
adds a pad to the 2-D space in percentage (by default set to 0.1) |
alpha |
a number between 0 and 1, to specify the transparency of the simulated line(s) |
resolution |
number of pixels the rasters are downsampled to (by default set to 500 pixels) |
A ggplot2 object.
plot2d(niclas)
plot2d(niclas)
Plot track(s) with a surface of a digital elevation model in three dimensions
plot3d( origTrack, simTrack = NULL, titleText = character(1), DEM = NULL, padding = 0.1, timesHeight = 10 )
plot3d( origTrack, simTrack = NULL, titleText = character(1), DEM = NULL, padding = 0.1, timesHeight = 10 )
origTrack |
a list containing data.frames with x,y,z coordinates or a data.frame |
simTrack |
a list containing data.frames with x,y,z coordinates or a data.frame |
titleText |
string with title of the plot |
DEM |
an object of type |
padding |
adds a pad to the 2-D space in percentage (by default set to 0.1) |
timesHeight |
multiply the height scale by a scalar (by default set to 10) |
Plots a plotly object
plot3d(niclas)
plot3d(niclas)
The function takes either one track or two tracks. The second track can be a list of tracks (eg. the output of n.sim.cond.3d), Then the densities of turn angle, lift angle and step length of all the simulations is taken. Additionally the autodifferences parameter can be set to true, then the densities of the autodifferences in turn angle, lift angle and step length are visualized.
plot3d.densities( track1, track2 = NULL, autodifferences = FALSE, scaleDensities = FALSE )
plot3d.densities( track1, track2 = NULL, autodifferences = FALSE, scaleDensities = FALSE )
track1 |
a list containing a data.frame with x,y,z coordinates or a data.frame |
track2 |
a list containing a data.frame with x,y,z coordinates or a data.frame |
autodifferences |
logical: should the densities of the autodifferences in turn angle, lift angle and step length are visualized. |
scaleDensities |
logical: should densities be scaled between 0 and 1, then sum of the area under the curve is not 1 anymore! |
A ggplot2 object.
plot3d.densities(niclas)
plot3d.densities(niclas)
If the layout is something like matrix(c(1,2,3,3), nrow=2, byrow=TRUE)
,
then plot 1 will go in the upper left, 2 will go in the upper right, and
3 will go all the way across the bottom.
plot3d.multiplot(..., plotlist = NULL, cols = 1, layout = NULL)
plot3d.multiplot(..., plotlist = NULL, cols = 1, layout = NULL)
... |
ggplot objects |
plotlist |
a list of ggplot objects |
cols |
number of columns in layout |
layout |
a matrix specifying the layout. If present, |
Nothing, plots the ggplot2 objects.
plot3d.multiplot(plot2d(niclas), plot2d(niclas), plot2d(niclas))
plot3d.multiplot(plot2d(niclas), plot2d(niclas), plot2d(niclas))
Creates a three dimensional scatterplot of the possibles next steps, based on the tldCube, which was extracted from a track.
plot3d.tldCube(tldCube)
plot3d.tldCube(tldCube)
tldCube |
tldCube; the ouptut from turnLiftStepHist or get.densities.3d |
Plots a plotly object
P <- get.track.densities.3d(niclas) suppressWarnings(plot3d.tldCube(P$tldCube))
P <- get.track.densities.3d(niclas) suppressWarnings(plot3d.tldCube(P$tldCube))
Plots a rasterLayer or rasterStack
plotRaster(r, title = character(0), centerColorBar = FALSE, ncol = NULL)
plotRaster(r, title = character(0), centerColorBar = FALSE, ncol = NULL)
r |
|
title |
title text of plot(s) |
centerColorBar |
logical: center colobar around 0 and use |
ncol |
number of columns to plot a stack, by default estimated by the square root |
Plots the rasters
plotRaster(dem)
plotRaster(dem)
Calculates the Q probability, representing the pull to the target. The number of steps on which the Q prob will be quantified is number of total segments less than one (the last step is defined by the target itself).
qProb.3d(sim, n.locs, parallel = FALSE, maxBin = 25)
qProb.3d(sim, n.locs, parallel = FALSE, maxBin = 25)
sim |
the result of sim.uncond.3d, or a data frame with at least x,y,z-coordinates, the arrival azimuth and the arrival gradient. |
n.locs |
number of total segments to be modeled, the length of the desired conditional empirical random walk |
parallel |
logical: run computations in parallel (n-1 cores)? Or numeric: the number of nodes (maximum: n - 1 cores) |
maxBin |
numeric scalar, maximum number of bins per dimension of the tld-cube (turnLiftStepHist) |
A list containing the Q - tldCubes for every step
qProb.3d(niclas, 3)
qProb.3d(niclas, 3)
Simulates n tracks with the geometrical properties of the original track, between the same start and end point.
reproduce.track.3d( track, n.sim = 1, parallel = FALSE, error = TRUE, DEM = NULL, BG = NULL, filterDeadEnds = TRUE, plot2d = FALSE, plot3d = FALSE, maxBin = 25, gradientDensity = TRUE )
reproduce.track.3d( track, n.sim = 1, parallel = FALSE, error = TRUE, DEM = NULL, BG = NULL, filterDeadEnds = TRUE, plot2d = FALSE, plot3d = FALSE, maxBin = 25, gradientDensity = TRUE )
track |
data.frame with x,y,z coordinates of the original track |
n.sim |
number of simulations that should be done |
parallel |
logical: run computations in parallel (n-1 cores)? Or numeric: the number of nodes (maximum: n - 1 cores) |
error |
logical: add error term to movement in simulation? |
DEM |
a raster containing a digital elevation model, covering the same extent as the track |
BG |
a raster influencing the probabilities. |
filterDeadEnds |
logical: Remove tracks that ended in a dead end? |
plot2d |
logical: plot tracks on 2-D plane? |
plot3d |
logical: plot tracks in 3-D? |
maxBin |
numeric scalar, maximum number of bins per dimension of the tld-cube (turnLiftStepHist) |
gradientDensity |
logical: Should a distribution of the gradient angle be extracted and used in the simulations (get.densities.3d)? |
A list or data.frame containing the simulated track(s) (CERW).
reproduce.track.3d(niclas[1:10, ])
reproduce.track.3d(niclas[1:10, ])
Exports a dataCube of type rasterStack
as Tiff image sequence.
Image sequences are a common structure to represent voxel data and
most of the specific software to visualize voxel data is able to read it (e.g. blender)
saveImageSlices(rStack, filename, dir, NaVal = 0)
saveImageSlices(rStack, filename, dir, NaVal = 0)
rStack |
rasterStack to be saved to Tiff image slices |
filename |
name of the image slices |
dir |
directory, where the slices should be stored |
NaVal |
numeric value that should represent NA values in the Tiff image, default is |
Saves the Tiff image files.
crws <- lapply(X = seq(1:100), FUN = function(X) { sim.crw.3d(nStep = 100, rTurn = 0.99, rLift = 0.99, meanStep = 0.1) }) points <- do.call("rbind", crws) extent <- raster::extent(c(-10, 10, -10, 10)) ud <- voxelCount(points, extent, xyRes = 5, zMin = -10, zMax = 10, standartize = TRUE ) saveImageSlices(ud, filename = "saveImageSlices_test", dir = tempdir())
crws <- lapply(X = seq(1:100), FUN = function(X) { sim.crw.3d(nStep = 100, rTurn = 0.99, rLift = 0.99, meanStep = 0.1) }) points <- do.call("rbind", crws) extent <- raster::extent(c(-10, 10, -10, 10)) ud <- voxelCount(points, extent, xyRes = 5, zMin = -10, zMax = 10, standartize = TRUE ) saveImageSlices(ud, filename = "saveImageSlices_test", dir = tempdir())
Converts a sf data.frame to a normal dataframe
sf2df.3d(track)
sf2df.3d(track)
track |
An object of type |
A data.frame.
sf2df.3d(track2sf.3d(niclas, CRS = 4326))
sf2df.3d(track2sf.3d(niclas, CRS = 4326))
Creates a conditional empirical random walk, with a specific starting and ending point, geometrically similar to the initial trajectory (extractMethod: raster overlay method can take "simple" or "bilinear")
sim.cond.3d( n.locs, start = c(0, 0, 0), end = start, a0, g0, densities, qProbs, error = FALSE, DEM = NULL, BG = NULL )
sim.cond.3d( n.locs, start = c(0, 0, 0), end = start, a0, g0, densities, qProbs, error = FALSE, DEM = NULL, BG = NULL )
n.locs |
length of the trajectory in locations |
start |
numeric vector of length 3 with the coordinates of the start point |
end |
numeric vector of length 3 with the coordinates of the end point |
a0 |
initial incoming heading in radian |
g0 |
initial incoming gradient/polar angle in radian |
densities |
list object returned by the get.densities.3d function |
qProbs |
list object returned by the qProb.3d function |
error |
logical: add random noise to the turn angle, lift angle and step length to account for errors measurements? |
DEM |
raster layer containing a digital elevation model, covering the area between start and end point |
BG |
a background raster layer that can be used to inform the choice of steps |
A trajectory in the form of data.frame
niclas <- track.properties.3d(niclas) n.locs <- 3 P <- get.track.densities.3d(niclas) f <- 1500 start <- Reduce(c, niclas[1, 1:3]) end <- Reduce(c, niclas[n.locs, 1:3]) a0 <- niclas$a[1] g0 <- niclas$g[1] uerw <- sim.uncond.3d( n.locs * f, start = start, a0 = a0, g0 = g0, densities = P ) Q <- qProb.3d(uerw, n.locs) sim.cond.3d( n.locs = n.locs, start = start, end = end, a0 = a0, g0 = g0, densities = P, qProbs = Q )
niclas <- track.properties.3d(niclas) n.locs <- 3 P <- get.track.densities.3d(niclas) f <- 1500 start <- Reduce(c, niclas[1, 1:3]) end <- Reduce(c, niclas[n.locs, 1:3]) a0 <- niclas$a[1] g0 <- niclas$g[1] uerw <- sim.uncond.3d( n.locs * f, start = start, a0 = a0, g0 = g0, densities = P ) Q <- qProb.3d(uerw, n.locs) sim.cond.3d( n.locs = n.locs, start = start, end = end, a0 = a0, g0 = g0, densities = P, qProbs = Q )
Simulation of a three dimensional Correlated Random Walk
sim.crw.3d(nStep, rTurn, rLift, meanStep, start = c(0, 0, 0))
sim.crw.3d(nStep, rTurn, rLift, meanStep, start = c(0, 0, 0))
nStep |
the number of steps of the simulated trajectory |
rTurn |
the correlation on the turn angle |
rLift |
the correlation of the lift angle |
meanStep |
the mean step length |
start |
a vector of length 3 containing the coordinates of the start point of the trajectory |
A trajectory in the form of data.frame
sim.crw.3d(nStep = 10, rTurn = 0.9, rLift = 0.9, meanStep = 1, start = c(0, 0, 0))
sim.crw.3d(nStep = 10, rTurn = 0.9, rLift = 0.9, meanStep = 1, start = c(0, 0, 0))
Creates a conditional empirical random walk in gliding mode, between a start and end point. The walk is performed on a MODE layer and, if provided, additionally on a background and digital elevation layer. The gliding is simulated with sim.cond.3d and soaring with sim.uncond.3d, therefore soaring is not restricted towards the target and can happen completly free as long as there are good thermal conditions. It is important to extract for every mode in the MODE raster layer a corresponding densities object with get.densities.3d and pass them to the function.
sim.glidingSoaring.3d( MODE, dGliding, dSoaring, qGliding, start = c(0, 0, 0), end = start, a0, g0, error = TRUE, smoothTransition = TRUE, glideRatio = 15, DEM = NULL, BG = NULL, verbose = FALSE )
sim.glidingSoaring.3d( MODE, dGliding, dSoaring, qGliding, start = c(0, 0, 0), end = start, a0, g0, error = TRUE, smoothTransition = TRUE, glideRatio = 15, DEM = NULL, BG = NULL, verbose = FALSE )
MODE |
raster layer containing the number/index of the mode, which should be used at each location |
dGliding |
density object returned by the get.densities.3d function for gliding mode |
dSoaring |
density object returned by the get.densities.3d function for soaring mode |
qGliding |
the Q probabilites for the steps in gliding mode (qProb.3d) |
start |
numeric vector of length 3 with the coordinates of the start point |
end |
numeric vector of length 3 with the coordinates of the end point |
a0 |
initial incoming heading in radian |
g0 |
initial incoming gradient/polar angle in radian |
error |
logical: add random noise to the turn angle, lift angle and step length to account for errors measurements? |
smoothTransition |
logical: should the transitions between soaring and the following gliding sections be smoothed? Recommended to avoid dead ends |
glideRatio |
ratio between vertical and horizontal movement, by default set to 15 meters forward movement per meter vertical movement |
DEM |
raster layer containing a digital elevation model, covering the area between start and end point |
BG |
a background raster layer that can be used to inform the choice of steps |
verbose |
logical: print current mode used? |
A 'soaring-gliding' trajectory in the form of data.frame
The MODE raster layer must be in the following structure: Gliding pixels have the value 1 and soaring pixel the values 2. NA
's are not allowed in the raster.
print("tbd.")
print("tbd.")
This function creates unconditional walks with prescribed empirical properties (turning angle, lift angle and step length and the auto-differences of them. It can be used for uncon- ditional walks or to seed the conditional walks with comparably long simulations. The conditional walk connecting a given start with a certain end point by a given number of steps needs an attraction term (the Q probability, see qProb.3d) to ensure that the target is approached and hit. In order to calculate the Q probability for each step the distribution of turns and lifts to target and the distribution of distance to target has to be known. They can be derived from the empirical data (ideally), or estimated from an unconditional process with the same properties. Creates a unconditional empirical random walk, with a specific starting point, geometrically similar to the initial trajectory.
sim.uncond.3d(n.locs, start = c(0, 0, 0), a0, g0, densities, error = TRUE)
sim.uncond.3d(n.locs, start = c(0, 0, 0), a0, g0, densities, error = TRUE)
n.locs |
the number of locations for the simulated track |
start |
vector indicating the start point |
a0 |
initial heading in radian |
g0 |
initial gradient/polar angle in radian |
densities |
list object returned by the get.densities.3d function |
error |
logical: add random noise to the turn angle, lift angle and step length to account for errors measurements? |
A 3 dimensional trajectory in the form of a data.frame
Simulations connecting start and end points
with more steps than 1/10th or more of the number of steps
of the empirical data should rather rely on simulated
unconditional walks with the same properties than on
the empirical data (factor = 1500
).
For a random initial heading a0 use:
sample(atan2(diff(coordinates(track)[,2]), diff(coordinates(track)[,1])),1)
sim.uncond.3d( 10, start = c(0, 0, 0), a0 = pi / 2, g0 = pi / 2, densities = get.track.densities.3d(niclas) )
sim.uncond.3d( 10, start = c(0, 0, 0), a0 = pi / 2, g0 = pi / 2, densities = get.track.densities.3d(niclas) )
The test simulates a CRW with given parameters and reconstructs it by using the eRTG3D
test.eRTG.3d( parallel = FALSE, returnResult = FALSE, plot2d = FALSE, plot3d = TRUE, plotDensities = TRUE )
test.eRTG.3d( parallel = FALSE, returnResult = FALSE, plot2d = FALSE, plot3d = TRUE, plotDensities = TRUE )
parallel |
logical: test running parallel? |
returnResult |
logical: return tracks generated? |
plot2d |
logical: plot tracks on 2-D plane? |
plot3d |
logical: plot tracks in 3-D? |
plotDensities |
logical: plot densities of turning angle, lift angle and step length? |
A list containing the original CRW and the simulated track (CERW).
test.eRTG.3d()
test.eRTG.3d()
Uses two-sample Kolmogorov-Smirnov test or the one-sample t-test to compare the geometric characteristics of the original track with the characteristics of the simulated track.
test.verification.3d(track1, track2, alpha = 0.05, plot = FALSE, test = "ks")
test.verification.3d(track1, track2, alpha = 0.05, plot = FALSE, test = "ks")
track1 |
data.frame or list of data.frames with x,y,z coordinates of the original track |
track2 |
data.frame or list of data.frames with x,y,z coordinates of the simulated track |
alpha |
scalar: significance level, default |
plot |
logical: plot the densities or differences of turn angle, lift angle and step length of the two tracks? |
test |
character: either |
Test objects of the 6 two-sample Kolmogorov-Smirnov test conducted.
By choosing test = "ttest"
a random sample, without replacement is taken from the longer track,
to shorten it to the length of the longer track. The order of the shorter track is also sampled randomly.
Then the two randomly ordered vectors of turn angles, lift angles and step lengths are substracted from each other.
If the both tracks stem from the same distributions the the mean deviatio should tend to towards zero, therefore the
difference is tested two-sided against mu = 0
with a one-sample t-test.
By setting test = "ks"
a two-sample Kolmogorov-Smirnov test is carried out on the distributions of turn angles,
lift angles and step lengths of the two tracks.
test.verification.3d(niclas, niclas)
test.verification.3d(niclas, niclas)
Extent of track(s)
track.extent(track, zAxis = FALSE)
track.extent(track, zAxis = FALSE)
track |
a list containing data.frames with x,y,z coordinates or a data.frame |
zAxis |
logical: return also the extent of the Z axis? |
Returns an extent object of the raster package in the 2–D case and a vector in the 3–D case.
track.extent(niclas, zAxis = TRUE)
track.extent(niclas, zAxis = TRUE)
Returns the properties (distances, azimuth, polar angle, turn angle & lift angle) of a track in three dimensions.
track.properties.3d(track)
track.properties.3d(track)
track |
data.frame with x,y,z coordinates |
The data.frame with track properties
track.properties.3d(niclas)
track.properties.3d(niclas)
The length of timeLag must be the the track's length minus 1 and represents the time passed between the fix point acquisition
track.split.3d(track, timeLag, lag = NULL, tolerance = NULL)
track.split.3d(track, timeLag, lag = NULL, tolerance = NULL)
track |
track data.frame with x, y and z coordinates |
timeLag |
a numeric vector with the time passed between the fix point acquisition |
lag |
NULL or a manually chosen lag |
tolerance |
NULL or a manually chosen tolerance |
A list containing the splitted tracks.
track.split.3d( niclas, timeLag = rep(1, nrow(niclas) - 1) + rnorm(nrow(niclas) - 1, mean = 0, sd = 0.25) )
track.split.3d( niclas, timeLag = rep(1, nrow(niclas) - 1) + rnorm(nrow(niclas) - 1, mean = 0, sd = 0.25) )
'sf, data.frame'
Converts a track to a 'sf, data.frame'
track2sf.3d(track, CRS = NA)
track2sf.3d(track, CRS = NA)
track |
eRTG3D track data.frame or a matrix |
CRS |
numeric, EPSG code of the CRS |
A track of type 'sf, data.frame'
.
track2sf.3d(niclas, 4326)
track2sf.3d(niclas, 4326)
Attention: Please use this function for CRS transformations,
since it is based on the st_transform from the sf package and therefore
supports CRS transformations in 3-D. Note: spTransform
from the sp
package
only supports transformations in the 2D plane, which will cause distortions
in the third dimension.
transformCRS.3d(track, fromCRS, toCRS)
transformCRS.3d(track, fromCRS, toCRS)
track |
data.frame with x,y,z coordinates |
fromCRS |
numeric, EPSG code of the current CRS |
toCRS |
numeric, EPSG code of the CRS to be converted in |
A data.frame containing x,y,z and variables.
transformCRS.3d(niclas, fromCRS = 2056, toCRS = 4326)
transformCRS.3d(niclas, fromCRS = 2056, toCRS = 4326)
Calculates the turn angle between every point in the track and the last point (target).
turn2target.3d(track)
turn2target.3d(track)
track |
a track data.frame containing x, y and z coordinates |
A numeric vector with the turn angles to target
turn2target.3d(niclas)
turn2target.3d(niclas)
Derives a three dimensional distribution of a turn angle, lift angle and step length, using the Freedman–Diaconis rule for estimating the number of bins.
turnLiftStepHist( turn, lift, step, printDims = TRUE, rm.zeros = TRUE, maxBin = 25 )
turnLiftStepHist( turn, lift, step, printDims = TRUE, rm.zeros = TRUE, maxBin = 25 )
turn |
numeric vector of turn angles |
lift |
numeric vector of lift angles |
step |
numeric vector of step lengths |
printDims |
logical: should dimensions of tld-Cube be messaged? |
rm.zeros |
logical: should combinations with zero probability be removed? |
maxBin |
numeric scalar, maximum number of bins per dimension of the tld-cube. |
A three dimensional histogram as data.frame
niclas <- track.properties.3d(niclas)[2:nrow(niclas), ] turnLiftStepHist(niclas$t, niclas$l, niclas$d)
niclas <- track.properties.3d(niclas)[2:nrow(niclas), ] turnLiftStepHist(niclas$t, niclas$l, niclas$d)
A rasterStack
object is created, representing the 3–D voxel cube.
The z axis is sliced into regular sections between the maximum and minimum value.
For every height slice a raster with points per cell counts is created. Additionally
the voxels can be standartized between 0 and 1.
voxelCount( points, extent, xyRes, zRes = xyRes, zMin, zMax, standartize = FALSE, verbose = FALSE )
voxelCount( points, extent, xyRes, zRes = xyRes, zMin, zMax, standartize = FALSE, verbose = FALSE )
points |
a x, y, z data.frame |
extent |
a raster extent object of the extent to create the rasters |
xyRes |
resolution in the ground plane of the created rasters |
zRes |
resolution in the z axis (by default |
zMin |
minimum z value |
zMax |
maximum height value |
standartize |
logical: standartize the values? |
verbose |
logical: print currently processed height band in raster stack? |
A rasterStack
object, representing the 3–D voxel cube.
voxelCount(niclas, raster::extent(dem), 100, 100, 1000, 1400, standartize = TRUE)
voxelCount(niclas, raster::extent(dem), 100, 100, 1000, 1400, standartize = TRUE)