Title: | Animal Movement Tools |
---|---|
Description: | Manage and analyze animal movement data. The functionality of 'amt' includes methods to calculate home ranges, track statistics (e.g. step lengths, speed, or turning angles), prepare data for fitting habitat selection analyses, and simulation of space-use from fitted step-selection functions. |
Authors: | Johannes Signer [aut, cre], Brian Smith [ctb], Bjoern Reineking [ctb], Ulrike Schlaegel [ctb], John Fieberg [ctb], Josh O'Brien [ctb], Bernardo Niebuhr [ctb], Alec Robitaille [ctb], Tal Avgar [ctb], Scott LaPoint [dtc] |
Maintainer: | Johannes Signer <[email protected]> |
License: | GPL-3 |
Version: | 0.3.0.0 |
Built: | 2024-11-13 15:26:27 UTC |
Source: | https://github.com/jmsigner/amt |
This file includes spatial data from 4 fisher (Pekania pennanti). These location data were collected via a 105g GPS tracking collar (manufactured by E-obs GmbH) and programmed to record the animal's location every 10 minutes, continuously. The data re projected in NAD84 (epsg: 5070). The data usage is permitted for exploratory purposes. For other purposes please get in contact (Scott LaPoint).
amt_fisher
amt_fisher
A tibble
with 14230 rows and 5 variables:
the x-coordinate
the y-coordinate
the timestamp
the sex of the animal
the id of the animal
the name of the animal
https://www.datarepository.movebank.org/handle/10255/move.330
For more information, contact Scott LaPoint [email protected]
A list with three entries that correspond to the following three layer: land use, elevation and population density.
amt_fisher_covar
amt_fisher_covar
A list with three where each entry is a SpatRast
.
https://lpdaac.usgs.gov/dataset_discovery/aster/aster_products_table
http://dup.esrin.esa.it/page_globcover.php
http://sedac.ciesin.columbia.edu/data/collection/gpw-v3/sets/browse
Exports a track to (multi)lines from the sf
package.
as_sf_lines(x, ...)
as_sf_lines(x, ...)
x |
|
... |
Further arguments, none implemented. |
A tibble
with a sfc
-column
Coerces a track to points from the sf
package.
as_sf_points(x, ...) ## S3 method for class 'steps_xy' as_sf_points(x, end = TRUE, ...)
as_sf_points(x, ...) ## S3 method for class 'steps_xy' as_sf_points(x, end = TRUE, ...)
x |
|
... |
Further arguments, none implemented. |
end |
|
A data data.frame
with a sfc
-column
Coerce other classes to a track_xy
.
as_track(x, ...) ## S3 method for class 'sfc_POINT' as_track(x, ...) ## S3 method for class 'steps_xyt' as_track(x, ...) ## S3 method for class 'data.frame' as_track(x, ...)
as_track(x, ...) ## S3 method for class 'sfc_POINT' as_track(x, ...) ## S3 method for class 'steps_xyt' as_track(x, ...) ## S3 method for class 'data.frame' as_track(x, ...)
x |
Object to be converted to a track. |
... |
Further arguments, none implemented. |
An object of class track_xy(t)
uhc_data
object to data.frame
Coerces uhc_data
from list
to data.frame
## S3 method for class 'uhc_data' as.data.frame(x, row.names = NULL, optional = FALSE, ...)
## S3 method for class 'uhc_data' as.data.frame(x, row.names = NULL, optional = FALSE, ...)
x |
|
row.names |
Included for consistency with generic
|
optional |
Included for consistency with generic
|
... |
Included for consistency with generic
|
This coercion aims to keep all of the information contained in
the uhc_data
list
in the resulting data.frame
representation. Factors
are converted to numeric, but the levels are retained in the column
"label"
.
Returns a data.frame
with columns:
var
: The name of the variable
x
: The x-coordinate of the density plot (the value of var
).
y
: The y-coordinate of the density plot (the probability density for
a numeric var
and the proportion for a factor var
).
dist
: The distribution represented. Either "U"
for used, "A"
for
available, or "S"
for sampled.
iter
: The iteration number if dist == "S"
.
label
: The label if var
is a factor.
Brian J. Smith
prep_uhc()
, conf_envelope()
Display available distributions for step lengths and turn angles.
available_distr(which_dist = "all", names_only = FALSE, ...)
available_distr(which_dist = "all", names_only = FALSE, ...)
which_dist |
|
names_only |
|
... |
none implemented. |
A tibble
with the purpose of the distribution (turn angles [ta] or step length [sl]) and the distribution name.
hr_kde_pi
wraps KernSmooth::dpik
to select bandwidth for kernel density estimation the plug-in-the-equation method in two dimensions.This function calculates bandwidths for kernel density estimation by wrapping KernSmooth::dpik
. If correct = TURE
, the bandwidth is transformed with power 5/6 to correct for using an univariate implementation for bivariate data (Gitzen et. al 2006).
hr_kde_pi(x, ...) ## S3 method for class 'track_xy' hr_kde_pi(x, rescale = "none", correct = TRUE, ...)
hr_kde_pi(x, ...) ## S3 method for class 'track_xy' hr_kde_pi(x, rescale = "none", correct = TRUE, ...)
x |
|
... |
Further arguments, none implemented. |
rescale |
|
correct |
Logical scalar that indicates whether or not the estimate should be correct for the two dimensional case. |
The bandwidth, the standardization method and correction.
Gitzen, R. A., Millspaugh, J. J., & Kernohan, B. J. (2006). Bandwidth selection for fixed-kernel analysis of animal utilization distributions. Journal of Wildlife Management, 70(5), 1334-1344.
KernSmooth::dpik
Calculate the reference bandwidth for kernel density home-range range estimates.
hr_kde_ref(x, ...) ## S3 method for class 'track_xy' hr_kde_ref(x, rescale = "none", ...)
hr_kde_ref(x, ...) ## S3 method for class 'track_xy' hr_kde_ref(x, rescale = "none", ...)
x |
|
... |
Further arguments, none implemented. |
rescale |
|
The estimated bandwidth in x and y direction.
Get bounding box of a track.
bbox(x, ...) ## S3 method for class 'track_xy' bbox(x, spatial = TRUE, buffer = NULL, ...) ## S3 method for class 'steps_xy' bbox(x, spatial = TRUE, buffer = NULL, ...)
bbox(x, ...) ## S3 method for class 'track_xy' bbox(x, spatial = TRUE, buffer = NULL, ...) ## S3 method for class 'steps_xy' bbox(x, spatial = TRUE, buffer = NULL, ...)
x |
|
... |
Further arguments, none implemented. |
spatial |
|
buffer |
|
If spatial = FALSE
a named vector of length four with the extent of the bounding box. Otherwise a SpatialPolygon
or a simple feature polygon with the bounding box.
data(deer) bbox(deer) bbox(deer, spatial = FALSE) bbox(deer, buffer = 100, spatial = FALSE) # For steps deer |> steps_by_burst() |> bbox(spatial = FALSE) deer |> steps_by_burst() |> bbox(buffer = 100, spatial = FALSE) deer |> steps_by_burst() |> random_steps() |> bbox(spatial = FALSE) # Further manipulations are possible deer |> bbox() |> sf::st_transform(4326)
data(deer) bbox(deer) bbox(deer, spatial = FALSE) bbox(deer, buffer = 100, spatial = FALSE) # For steps deer |> steps_by_burst() |> bbox(spatial = FALSE) deer |> steps_by_burst() |> bbox(buffer = 100, spatial = FALSE) deer |> steps_by_burst() |> random_steps() |> bbox(spatial = FALSE) # Further manipulations are possible deer |> bbox() |> sf::st_transform(4326)
w(x)
Calculates the value of the exponential habitat selection function
calc_w(f, b, newdata)
calc_w(f, b, newdata)
f |
|
b |
|
newdata |
|
This is actually like to be w(x) * phi(x) for an iSSF.
Calculates squared displacement rate for a given speed and duration
calculate_sdr(speed = 50, time, speed_unit = c("km/h", "m/s"))
calculate_sdr(speed = 50, time, speed_unit = c("km/h", "m/s"))
speed |
|
time |
|
speed_unit |
|
Returns a numeric vector (of length 1) with the SDR in m^2/s
.
Johannes Signer and Brian J. Smith
# Assume a cheetah can sprint 100 km/h for 60 seconds calculate_sdr(speed = 100, time = seconds(60), speed_unit = "km/h") # 46296.3 m^2/s # What is expected displacement in 1 h at that SDR? get_displacement(46296.3, hours(1)) # 12909.95 m = 12.9 km/h (much slower than sprint speed!)
# Assume a cheetah can sprint 100 km/h for 60 seconds calculate_sdr(speed = 100, time = seconds(60), speed_unit = "km/h") # 46296.3 m^2/s # What is expected displacement in 1 h at that SDR? get_displacement(46296.3, hours(1)) # 12909.95 m = 12.9 km/h (much slower than sprint speed!)
Calculate the centroid of a track.
centroid(x, ...) ## S3 method for class 'track_xy' centroid(x, spatial = FALSE, ...)
centroid(x, ...) ## S3 method for class 'track_xy' centroid(x, spatial = FALSE, ...)
x |
|
... |
Further arguments, none implemented. |
spatial |
|
The centroid of a track as numeric vector if spatial = FALSE
, otherwise as SpatialPoints
.
data(deer) centroid(deer)
data(deer) centroid(deer)
Calculates speed
check_time_unit(tu)
check_time_unit(tu)
tu |
The |
Calculate Change in NSD
Calculates change in NSD
Check time_unit
Parameter
Internal function to check time_unit
parameter in various cleaning functions.
Several other packages provides methods to analyze movement data, and amt
provides coercion methods to some packages.
as_sf(x, ...) ## S3 method for class 'steps_xy' as_sf(x, end = TRUE, ...) as_sp(x, ...) as_ltraj(x, ...) ## S3 method for class 'track_xy' as_ltraj(x, id = "animal_1", ...) ## S3 method for class 'track_xyt' as_ltraj(x, ...) as_telemetry(x, ...) ## S3 method for class 'track_xyt' as_telemetry(x, ...) as_moveHMM(x, ...) ## S3 method for class 'track_xy' as_moveHMM(x, ...)
as_sf(x, ...) ## S3 method for class 'steps_xy' as_sf(x, end = TRUE, ...) as_sp(x, ...) as_ltraj(x, ...) ## S3 method for class 'track_xy' as_ltraj(x, id = "animal_1", ...) ## S3 method for class 'track_xyt' as_ltraj(x, ...) as_telemetry(x, ...) ## S3 method for class 'track_xyt' as_telemetry(x, ...) as_moveHMM(x, ...) ## S3 method for class 'track_xy' as_moveHMM(x, ...)
x |
|
... |
Further arguments, none implemented. |
end |
|
id |
|
An object of the class to which coercion is performed to.
uhc_data_frame
Simplifies sampled distributions in a uhc_data_frame
to confidence envelopes
conf_envelope(x, levels = c(0.95, 1))
conf_envelope(x, levels = c(0.95, 1))
x |
|
levels |
|
This can dramatically improve plotting time for UHC plots by simplifying the many sampled lines down to the boundaries of a polygon.
Returns a data.frame
with columns:
var
: The name of the variable
x
: The x-coordinate of the density plot (the value of var
).
label
: If var
is a factor
, the label for the value given by x
.
U
: The y-coordinate of the density plot for the use distribution.
A
: The y-coordinate of the density plot for the availability distribution.
CI*_lwr
: The lower bound of the confidence envelope for the corresponding
confidence level.
CI*_upr
: The upper bound of the confidence envelope for the corresponding
confidence level.
Brian J. Smith
prep_uhc()
, plot.uhc_envelopes()
Converts angles to radians
as_rad(x) as_degree(x)
as_rad(x) as_degree(x)
x |
|
A numeric vector with the converted angles.
as_rad(seq(-180, 180, 30)) # The default unit of turning angles is rad. data(deer) deer |> steps() |> mutate(ta_ = as_degree(ta_))
as_rad(seq(-180, 180, 30)) # The default unit of turning angles is rad. data(deer) deer |> steps() |> mutate(ta_ = as_degree(ta_))
Coordinates of a track.
coords(x, ...)
coords(x, ...)
x |
|
... |
Further arguments, none implemented. |
[tibble]
The coordinates.
data(deer) coords(deer)
data(deer) coords(deer)
Calculate the cumulative utilization distribution (UD).
hr_cud(x, ...) ## S3 method for class 'SpatRaster' hr_cud(x, ...)
hr_cud(x, ...) ## S3 method for class 'SpatRaster' hr_cud(x, ...)
x |
|
... |
Further arguments, none implemented. |
[RasterLayer]
The cumulative UD.
This function is typically used to obtain isopleths.
826 GPS relocations of one red deer from northern Germany. The data is
already resampled to a regular time interval of 6 hours and the coordinate
reference system is transformed to epsg:3035
.
deer
deer
A track_xyt
the x-coordinate
the y-coordinate
the timestamp
the burst a particular points belongs to.
Verein für Wildtierforschung Dresden und Göttingen e.V.
Difference in x and y coordinates.
diff_x(x, ...) diff_y(x, ...)
diff_x(x, ...) diff_y(x, ...)
x |
|
... |
Further arguments, none implemented. |
Numeric vector
Name of step-length distribution and turn-angle distribution
sl_distr_name(x, ...) ## S3 method for class 'random_steps' sl_distr_name(x, ...) ## S3 method for class 'fit_clogit' sl_distr_name(x, ...) ta_distr_name(x, ...) ta_distr_name(x, ...) ## S3 method for class 'random_steps' ta_distr_name(x, ...) ## S3 method for class 'fit_clogit' ta_distr_name(x, ...)
sl_distr_name(x, ...) ## S3 method for class 'random_steps' sl_distr_name(x, ...) ## S3 method for class 'fit_clogit' sl_distr_name(x, ...) ta_distr_name(x, ...) ta_distr_name(x, ...) ## S3 method for class 'random_steps' ta_distr_name(x, ...) ## S3 method for class 'fit_clogit' ta_distr_name(x, ...)
x |
Random steps or fitted model |
... |
None implemented. |
Character vector of length 1.
make_distributions
creates a distribution suitable for using it with integrated step selection functions
make_distribution(name, params, vcov = NULL, ...) make_exp_distr(rate = 1) make_hnorm_distr(sd = 1) make_lnorm_distr(meanlog = 0, sdlog = 1) make_unif_distr(min = -pi, max = pi) make_vonmises_distr(kappa = 1, vcov = NULL) make_gamma_distr(shape = 1, scale = 1, vcov = NULL)
make_distribution(name, params, vcov = NULL, ...) make_exp_distr(rate = 1) make_hnorm_distr(sd = 1) make_lnorm_distr(meanlog = 0, sdlog = 1) make_unif_distr(min = -pi, max = pi) make_vonmises_distr(kappa = 1, vcov = NULL) make_gamma_distr(shape = 1, scale = 1, vcov = NULL)
name |
|
params |
|
vcov |
|
... |
none implemented. |
rate |
|
sd |
|
meanlog |
|
sdlog |
|
min |
|
max |
|
kappa |
|
shape , scale
|
|
A list of class amt_distr
that contains the name (name
) and parameters (params
) of a distribution.
Obtain the extent of a track in x
y
or both directions
extent_x(x, ...) extent_y(x, ...) extent_both(x, ...) extent_max(x, ...)
extent_x(x, ...) extent_y(x, ...) extent_both(x, ...) extent_max(x, ...)
x |
|
... |
Further arguments, none implemented. |
Numeric vector with the extent.
Extract the covariate values at relocations, or at the beginning or end of steps.
extract_covariates(x, ...) ## S3 method for class 'track_xy' extract_covariates(x, covariates, ...) ## S3 method for class 'random_points' extract_covariates(x, covariates, ...) ## S3 method for class 'steps_xy' extract_covariates(x, covariates, where = "end", ...) extract_covariates_along(x, ...) ## S3 method for class 'steps_xy' extract_covariates_along(x, covariates, ...) extract_covariates_var_time(x, ...) ## S3 method for class 'track_xyt' extract_covariates_var_time( x, covariates, when = "any", max_time, name_covar = "time_var_covar", ... ) ## S3 method for class 'steps_xyt' extract_covariates_var_time( x, covariates, when = "any", max_time, name_covar = "time_var_covar", where = "end", ... )
extract_covariates(x, ...) ## S3 method for class 'track_xy' extract_covariates(x, covariates, ...) ## S3 method for class 'random_points' extract_covariates(x, covariates, ...) ## S3 method for class 'steps_xy' extract_covariates(x, covariates, where = "end", ...) extract_covariates_along(x, ...) ## S3 method for class 'steps_xy' extract_covariates_along(x, covariates, ...) extract_covariates_var_time(x, ...) ## S3 method for class 'track_xyt' extract_covariates_var_time( x, covariates, when = "any", max_time, name_covar = "time_var_covar", ... ) ## S3 method for class 'steps_xyt' extract_covariates_var_time( x, covariates, when = "any", max_time, name_covar = "time_var_covar", where = "end", ... )
x |
|
... |
Additional arguments passed to |
covariates |
|
where |
|
when |
|
max_time |
|
name_covar |
|
extract_covariates_along
extracts the covariates along a straight line between the start and the end point of a (random) step. It returns a list, which in most cases will have to be processed further.
A tibble
with additional columns for covariate values.
data(deer) sh_forest <- get_sh_forest() mini_deer <- deer[1:20, ] mini_deer |> extract_covariates(sh_forest) mini_deer |> steps() |> extract_covariates(sh_forest) # Illustration of extracting covariates along the a step mini_deer |> steps() |> random_steps() |> extract_covariates(sh_forest) |> # extract at the endpoint (\(.) mutate(., for_path = extract_covariates_along(., sh_forest)))() |> # 1 = forest, lets calc the fraction of forest along the path mutate(for_per = purrr::map_dbl(for_path, function(x) mean(x == 1)))
data(deer) sh_forest <- get_sh_forest() mini_deer <- deer[1:20, ] mini_deer |> extract_covariates(sh_forest) mini_deer |> steps() |> extract_covariates(sh_forest) # Illustration of extracting covariates along the a step mini_deer |> steps() |> random_steps() |> extract_covariates(sh_forest) |> # extract at the endpoint (\(.) mutate(., for_path = extract_covariates_along(., sh_forest)))() |> # 1 = forest, lets calc the fraction of forest along the path mutate(for_per = purrr::map_dbl(for_path, function(x) mean(x == 1)))
uhc_data
objectSubset a uhc_data
object
## S3 method for class 'uhc_data' x[i]
## S3 method for class 'uhc_data' x[i]
x |
|
i |
|
Only retain bursts with a minimum number (= min_n
) of relocations.
filter_min_n_burst(x, ...) ## S3 method for class 'track_xy' filter_min_n_burst(x, min_n = 3, ...)
filter_min_n_burst(x, ...) ## S3 method for class 'track_xy' filter_min_n_burst(x, min_n = 3, ...)
x |
|
... |
Further arguments, none implemented. |
min_n |
|
A tibble
of class track_xy(t)
.
This function is a wrapper around survival::clogit
, making it usable in a piped workflow.
fit_clogit(data, formula, more = NULL, summary_only = FALSE, ...) fit_ssf(data, formula, more = NULL, summary_only = FALSE, ...) fit_issf(data, formula, more = NULL, summary_only = FALSE, ...)
fit_clogit(data, formula, more = NULL, summary_only = FALSE, ...) fit_ssf(data, formula, more = NULL, summary_only = FALSE, ...) fit_issf(data, formula, more = NULL, summary_only = FALSE, ...)
data |
|
formula |
|
more |
|
summary_only |
|
... |
Additional arguments, passed to |
A list with the following entries
model: The model output.
sl_: The step length distribution.
ta_: The turn angle distribution.
ctmm
Fit a continuous time movement model with ctmm
fit_ctmm(x, model, uere = NULL, ...)
fit_ctmm(x, model, uere = NULL, ...)
x |
|
model |
|
uere |
User Equivalent Range Error, see |
... |
Additional parameters passed to |
An object of class ctmm
from the package ctmm.
C. H. Fleming, J. M. Calabrese, T. Mueller, K.A. Olson, P. Leimgruber, W. F. Fagan, “From fine-scale foraging to home ranges: A semi-variance approach to identifying movement modes across spatiotemporal scales”, The American Naturalist, 183:5, E154-E167 (2014).
data(deer) mini_deer <- deer[1:20, ] m1 <- fit_ctmm(mini_deer, "iid") summary(m1)
data(deer) mini_deer <- deer[1:20, ] m1 <- fit_ctmm(mini_deer, "iid") summary(m1)
Wrapper to fit a distribution to data. Currently implemented distributions
are the exponential distribution (exp
), the gamma distribution (gamma
)
and the von Mises distribution (vonmises
).
fit_distr(x, dist_name, na.rm = TRUE)
fit_distr(x, dist_name, na.rm = TRUE)
x |
|
dist_name |
|
na.rm |
|
An amt_distr
object, which consists of a list with the name
of
the distribution and its parameters (saved in params
).
set.seed(123) dat <- rexp(1e3, 2) fit_distr(dat, "exp")
set.seed(123) dat <- rexp(1e3, 2) fit_distr(dat, "exp")
This function is a wrapper around stats::glm
for a piped workflows.
fit_logit(data, formula, ...) fit_rsf(data, formula, ...)
fit_logit(data, formula, ...) fit_rsf(data, formula, ...)
data |
|
formula |
|
... |
Further arguments passed to |
A list with the model output.
Flags defunct clusters at the end of a track
flag_defunct_clusters(x, zeta, eta, theta, ...) ## S3 method for class 'track_xyt' flag_defunct_clusters(x, zeta, eta, theta, ...)
flag_defunct_clusters(x, zeta, eta, theta, ...) ## S3 method for class 'track_xyt' flag_defunct_clusters(x, zeta, eta, theta, ...)
x |
|
zeta |
|
eta |
|
theta |
|
... |
Addtional arguments. None currently implemented. |
Locations at the end of a trajectory may represent a dropped collar or an animal mortality. In some cases, the device may be recording locations for quite some time that are not biologically meaningful. This function aims to flag those locations at the end of the trajectory that belong to a mortality (or similar) cluster. The first location at the cluster remains unflagged, but all subsequent locations are flagged.
The algorithm detects steps that represent zero movement, within a precision
threshold given by zeta
. That is, if zeta = 5
(units determined by CRS;
typically meters), all points that differ by less than 5 will be considered
zero movement. Consecutive steps of zero movement (within the tolerance) form
a cluster. The parameter eta
gives the cutoff for the minimum number of
zero steps to be considered a cluster. Finally, the algorithm requires that
clusters persist without a non-zero step for a minimum amount of time, given
by theta
.
Returns x
(a track_xyt
) with a flagging column added
(x$defunct_cluster_
).
Brian J. Smith and Johannes Signer, based on code by Tal Avgar
flag_duplicates()
,
flag_fast_steps()
,
flag_roundtrips()
Flags locations with duplicate timestamps by DOP and distance
flag_duplicates(x, gamma, time_unit = "mins", DOP = "dop", ...) ## S3 method for class 'track_xyt' flag_duplicates(x, gamma, time_unit = "mins", DOP = "dop", ...)
flag_duplicates(x, gamma, time_unit = "mins", DOP = "dop", ...) ## S3 method for class 'track_xyt' flag_duplicates(x, gamma, time_unit = "mins", DOP = "dop", ...)
x |
|
gamma |
|
time_unit |
|
DOP |
|
... |
Additional arguments. None currently implemented. |
Locations are considered duplicates if their timestamps are within
gamma
of each other. However, the function runs sequentially through the
track object, so that only timestamps after the focal point are flagged as
duplicates (and thus removed from further consideration). E.g., if
gamma = minutes(5)
, then all locations with timestamp within 5 minutes
after the focal location will be considered duplicates.
When duplicates are found, (1) the location with the lowest dilution of precision
(given by DOP
column) is kept. If there are multiple duplicates with equally
low DOP, then (2) the one closest in space to previous location is kept. In
the event of exact ties in DOP and distance, (3) the first location is kept.
This is unlikely unless there are exact coordinate duplicates.
In the case that the first location in a trajectory has a duplicate, there is no previous location with which to calculate a distance. In that case, the algorithm skips to (3) and keeps the first location.
In the event your data.frame
does not have a DOP column, you can insert a dummy
with constant values such that all duplicates will tie, and distance will be
the only criterion (e.g., x$dop <- 1
). In the event you do have an alternate
measure of precision where larger numbers are more precise (e.g., number of
satellites), simply multiply that metric by -1
and pass it as if it were DOP.
Internally, the function drops duplicates as it works sequentially through the
data.frame
. E.g., if location 5 was considered a duplicate of location 4 –
and location 4 was higher quality – then location 5 would be dropped. The
function would then move on to location 6 (since 5 was already dropped).
However, the object returned to the user has all the original rows of x
–
i.e., locations are flagged rather than removed.
Returns x
(a track_xyt
) with a flagging column added (x$duplicate_
).
Brian J. Smith, based on code by Johannes Signer and Tal Avgar
flag_fast_steps()
,
flag_roundtrips()
,
flag_defunct_clusters()
Flags locations that imply SDR exceeding a threshold
flag_fast_steps(x, delta, time_unit = "secs", ...) ## S3 method for class 'track_xyt' flag_fast_steps(x, delta, time_unit = "secs", ...)
flag_fast_steps(x, delta, time_unit = "secs", ...) ## S3 method for class 'track_xyt' flag_fast_steps(x, delta, time_unit = "secs", ...)
x |
|
delta |
|
time_unit |
|
... |
Addtional arguments. None currently implemented. |
Locations are flagged if the SDR from the previous location to the
current location exceeds delta
. Internally, flagged locations are dropped
from future consideration.
The time_unit
should be the same time unit with which the SDR threshold
was calculated. SDR is typically calculated in m^2/s
, so time_unit
defaults
to "secs"
. The spatial unit is determined by the CRS, which should typically
be in meters.
Returns x
(a track_xyt
) with a flagging column added
(x$fast_step_
).
Brian J. Smith, based on code by Johannes Signer and Tal Avgar
flag_duplicates()
, flag_roundtrips()
,
flag_defunct_clusters()
Flags locations that imply fast round trips
flag_roundtrips(x, delta, epsilon, time_unit = "secs", ...) ## S3 method for class 'track_xyt' flag_roundtrips(x, delta, epsilon, time_unit = "secs", ...)
flag_roundtrips(x, delta, epsilon, time_unit = "secs", ...) ## S3 method for class 'track_xyt' flag_roundtrips(x, delta, epsilon, time_unit = "secs", ...)
x |
|
delta |
|
epsilon |
|
time_unit |
|
... |
Addtional arguments. None currently implemented. |
Locations implying a single fast step can be flagged using
flag_fast_steps()
. However, it is more likely that a single
location is imprecise if it implies an unrealistically fast out-and-back round
trip. In that case, the user might be willing to scale the threshold SDR.
In this function, delta
gives the base SDR and epsilon
is the scaling
factor, such that locations are considered for flagging if the SDR from the
previous location (location i - 1) to the focal location (i) [call it sdr1
]
and the focal location (i) to the next location (i + 1) [call it sdr2
] both
have SDR > delta/epsilon
.
In that case, the SDR from the previous location (i - 1) to the next location
(i + 1) is computed; i.e., the SDR assuming we omit the focal location (i)
[call it sdr3
]. The remaining locations should be closer together than
they are to the omitted location. Thus the focal location is flagged if
(sdr1 > epsilon * sdr3) & (sdr2 > epsilon * sdr3)
.
Note that epsilon
both decreases delta
in the out-and-back case and
increases sdr3
(between the remaining neighbors).
Internally, flagged locations are dropped from future consideration.
The time_unit
should be the same time unit with which the SDR threshold
was calculated. SDR is typically calculated in m^2/s
, so time_unit
defaults
to "secs"
. The spatial unit is determined by the CRS, which should typically
be in meters. The epsilon
parameter is unitless.
Returns x
(a track_xyt
) with a flagging column added
(x$fast_roundtrip_
).
Brian J. Smith, based on code by Johannes Signer and Tal Avgar
flag_duplicates()
,
flag_fast_steps()
,
flag_defunct_clusters()
Function that returns the start (from
), end (to
), and the duration (from_to
) of a track.
from_to(x, ...) ## S3 method for class 'track_xyt' from_to(x, ...) from(x, ...) ## S3 method for class 'track_xyt' from(x, ...) to(x, ...) ## S3 method for class 'track_xyt' to(x, ...)
from_to(x, ...) ## S3 method for class 'track_xyt' from_to(x, ...) from(x, ...) ## S3 method for class 'track_xyt' from(x, ...) to(x, ...) ## S3 method for class 'track_xyt' to(x, ...)
x |
|
... |
Further arguments, none implemented. |
A vector of class POSIXct
.
The current version of terra
(1.7.12) requires SpatRast
ers to be wrapped in order to be saved locally. This function unwraps the covariates for the fisher data and returns a list.
get_amt_fisher_covars()
get_amt_fisher_covars()
A list with covariates
Returns the proj4string
of an object.
get_crs(x, ...)
get_crs(x, ...)
x |
|
... |
Further arguments, none implemented. |
The proj4string
of the CRS.
data(deer) get_crs(deer)
data(deer) get_crs(deer)
Calculates expected displacement for a given SDR and time span
get_displacement(delta, time_span)
get_displacement(delta, time_span)
delta |
|
time_span |
|
Returns a numeric vector (of length 1) with the expected displacement in meters.
Johannes Signer and Brian J. Smith
Obtain the step length and/or turn angle distributions from random steps or a fitted model.
sl_distr(x, ...) ## S3 method for class 'random_steps' sl_distr(x, ...) ## S3 method for class 'fit_clogit' sl_distr(x, ...) ta_distr(x, ...) ## S3 method for class 'random_steps' ta_distr(x, ...) ## S3 method for class 'fit_clogit' ta_distr(x, ...)
sl_distr(x, ...) ## S3 method for class 'random_steps' sl_distr(x, ...) ## S3 method for class 'fit_clogit' sl_distr(x, ...) ta_distr(x, ...) ## S3 method for class 'random_steps' ta_distr(x, ...) ## S3 method for class 'fit_clogit' ta_distr(x, ...)
x |
Random steps or fitted model |
... |
None implemented. |
An amt distribution
Helper function to get the maximum distance from a fitted model.
get_max_dist(x, ...) ## S3 method for class 'fit_clogit' get_max_dist(x, p = 0.99, ...)
get_max_dist(x, ...) ## S3 method for class 'fit_clogit' get_max_dist(x, p = 0.99, ...)
x |
|
... |
Further arguments, none implemented. |
p |
|
The current version of terra
(1.7.12) requires SpatRast
ers to be wrapped in order to be saved locally. This function unwraps the the forest layer and returns a SpatRast
.
get_sh_forest()
get_sh_forest()
A SpatRast
with forest cover.
Checks if an object has a CRS.
has_crs(x, ...)
has_crs(x, ...)
x |
|
... |
Further arguments, none implemented. |
Logic vector of length 1.
data(deer) has_crs(deer)
data(deer) has_crs(deer)
Functions to calculate animal home ranges from a track_xy*
. hr_mcp
, hr_kde
, and hr_locoh
calculate the minimum convex
polygon, kernel density, and local convex hull home range respectively.
hr_akde(x, ...) ## S3 method for class 'track_xyt' hr_akde( x, model = fit_ctmm(x, "iid"), keep.data = TRUE, trast = make_trast(x), levels = 0.95, wrap = FALSE, ... ) hr_kde(x, ...) ## S3 method for class 'track_xy' hr_kde( x, h = hr_kde_ref(x), trast = make_trast(x), levels = 0.95, keep.data = TRUE, wrap = FALSE, ... ) hr_locoh(x, ...) ## S3 method for class 'track_xy' hr_locoh( x, n = 10, type = "k", levels = 0.95, keep.data = TRUE, rand_buffer = 1e-05, ... ) hr_mcp(x, ...) hr_od(x, ...)
hr_akde(x, ...) ## S3 method for class 'track_xyt' hr_akde( x, model = fit_ctmm(x, "iid"), keep.data = TRUE, trast = make_trast(x), levels = 0.95, wrap = FALSE, ... ) hr_kde(x, ...) ## S3 method for class 'track_xy' hr_kde( x, h = hr_kde_ref(x), trast = make_trast(x), levels = 0.95, keep.data = TRUE, wrap = FALSE, ... ) hr_locoh(x, ...) ## S3 method for class 'track_xy' hr_locoh( x, n = 10, type = "k", levels = 0.95, keep.data = TRUE, rand_buffer = 1e-05, ... ) hr_mcp(x, ...) hr_od(x, ...)
x |
|
... |
Further arguments, none implemented. |
model |
A continuous time movement model. This can be fitted either with |
keep.data |
|
trast |
|
levels |
|
wrap |
|
h |
|
n |
|
type |
|
rand_buffer |
|
A hr
-estimate.
Worton, B. J. (1989). Kernel methods for estimating the utilization distribution in home-range studies. Ecology, 70(1), 164-168. C. H. Fleming, W. F. Fagan, T. Mueller, K. A. Olson, P. Leimgruber, J. M. Calabrese, “Rigorous home-range estimation with movement data: A new autocorrelated kernel-density estimator”, Ecology, 96:5, 1182-1188 (2015).
Fleming, C. H., Fagan, W. F., Mueller, T., Olson, K. A., Leimgruber, P., & Calabrese, J. M. (2016). Estimating where and how animals travel: an optimal framework for path reconstruction from autocorrelated tracking data. Ecology, 97(3), 576-582.
data(deer) mini_deer <- deer[1:100, ] # MCP --------------------------------------------------------------------- mcp1 <- hr_mcp(mini_deer) hr_area(mcp1) # calculated MCP at different levels mcp1 <- hr_mcp(mini_deer, levels = seq(0.3, 1, 0.1)) hr_area(mcp1) # CRS are inherited get_crs(mini_deer) mcps <- hr_mcp(mini_deer, levels = c(0.5, 0.95, 1)) has_crs(mcps) # Kernel density estimaiton (KDE) ----------------------------------------- kde1 <- hr_kde(mini_deer) hr_area(kde1) get_crs(kde1) # akde data(deer) mini_deer <- deer[1:20, ] ud1 <- hr_akde(mini_deer) # uses an iid ctmm ud2 <- hr_akde(mini_deer, model = fit_ctmm(deer, "ou")) # uses an OU ctmm # od data(deer) ud1 <- hr_od(deer) # uses an iid ctmm ud2 <- hr_akde(deer, model = fit_ctmm(deer, "ou")) # uses an OU ctmm
data(deer) mini_deer <- deer[1:100, ] # MCP --------------------------------------------------------------------- mcp1 <- hr_mcp(mini_deer) hr_area(mcp1) # calculated MCP at different levels mcp1 <- hr_mcp(mini_deer, levels = seq(0.3, 1, 0.1)) hr_area(mcp1) # CRS are inherited get_crs(mini_deer) mcps <- hr_mcp(mini_deer, levels = c(0.5, 0.95, 1)) has_crs(mcps) # Kernel density estimaiton (KDE) ----------------------------------------- kde1 <- hr_kde(mini_deer) hr_area(kde1) get_crs(kde1) # akde data(deer) mini_deer <- deer[1:20, ] ud1 <- hr_akde(mini_deer) # uses an iid ctmm ud2 <- hr_akde(mini_deer, model = fit_ctmm(deer, "ou")) # uses an OU ctmm # od data(deer) ud1 <- hr_od(deer) # uses an iid ctmm ud2 <- hr_akde(deer, model = fit_ctmm(deer, "ou")) # uses an OU ctmm
Obtain the area of a home-range estimate, possible at different isopleth levels.
hr_area(x, ...) ## S3 method for class 'hr' hr_area(x, units = FALSE, ...) ## S3 method for class 'SpatRaster' hr_area(x, level = 0.95, ...) ## S3 method for class 'akde' hr_area(x, units = FALSE, ...)
hr_area(x, ...) ## S3 method for class 'hr' hr_area(x, units = FALSE, ...) ## S3 method for class 'SpatRaster' hr_area(x, level = 0.95, ...) ## S3 method for class 'akde' hr_area(x, units = FALSE, ...)
x |
An object of class |
... |
Further arguments, none implemented. |
units |
|
level |
The level at which the area will be calculated. |
A tibble
with the home-range level and the area.
Obtain the isopleths of a home-range estimate, possible at different isopleth levels.
hr_isopleths(x, ...) ## S3 method for class 'PackedSpatRaster' hr_isopleths(x, levels, descending = TRUE, ...) ## S3 method for class 'SpatRaster' hr_isopleths(x, levels, descending = TRUE, ...) ## S3 method for class 'mcp' hr_isopleths(x, descending = TRUE, ...) ## S3 method for class 'locoh' hr_isopleths(x, descending = TRUE, ...) ## S3 method for class 'hr_prob' hr_isopleths(x, descending = TRUE, ...) ## S3 method for class 'akde' hr_isopleths(x, conf.level = 0.95, descending = TRUE, ...)
hr_isopleths(x, ...) ## S3 method for class 'PackedSpatRaster' hr_isopleths(x, levels, descending = TRUE, ...) ## S3 method for class 'SpatRaster' hr_isopleths(x, levels, descending = TRUE, ...) ## S3 method for class 'mcp' hr_isopleths(x, descending = TRUE, ...) ## S3 method for class 'locoh' hr_isopleths(x, descending = TRUE, ...) ## S3 method for class 'hr_prob' hr_isopleths(x, descending = TRUE, ...) ## S3 method for class 'akde' hr_isopleths(x, conf.level = 0.95, descending = TRUE, ...)
x |
An object of class |
... |
Further arguments, none implemented. |
levels |
|
descending |
|
conf.level |
The confidence level for isopleths for |
A tibble
with the home-range level and a simple feature columns with the isoploth as multipolygon.
Use least square cross validation (lscv) to estimate bandwidth for kernel home-range estimation.
hr_kde_lscv( x, range = do.call(seq, as.list(c(hr_kde_ref(x) * c(0.1, 2), length.out = 100))), which_min = "global", rescale = "none", trast = make_trast(x) )
hr_kde_lscv( x, range = do.call(seq, as.list(c(hr_kde_ref(x) * c(0.1, 2), length.out = 100))), which_min = "global", rescale = "none", trast = make_trast(x) )
x |
|
range |
numeric vector with different candidate h values. |
which_min |
A character indicating if the |
rescale |
|
trast |
A template raster. |
hr_kde_lscv
calculates least square cross validation bandwidth. This implementation is based on Seaman and Powell (1996). If whichMin
is "global"
the global minimum is returned, else the local minimum with the largest candidate bandwidth is returned.
vector
of length two.
Seaman, D. E., & Powell, R. A. (1996). An evaluation of the accuracy of kernel density estimators for home range analysis. Ecology, 77(7), 2075-2085.
Use two dimensional reference bandwidth to select a bandwidth for kernel density estimation. Find the smallest value for bandwidth (h) that results in n polygons (usually n=1) contiguous polygons at a given level.
hr_kde_ref_scaled( x, range = hr_kde_ref(x)[1] * c(0.01, 1), trast = make_trast(x), num.of.parts = 1, levels = 0.95, tol = 0.1, max.it = 500L )
hr_kde_ref_scaled( x, range = hr_kde_ref(x)[1] * c(0.01, 1), trast = make_trast(x), num.of.parts = 1, levels = 0.95, tol = 0.1, max.it = 500L )
x |
A |
range |
Numeric vector, indicating the lower and upper bound of the search range. If |
trast |
A template |
num.of.parts |
Numeric numeric scalar, indicating the number of contiguous polygons desired. This will usually be one. |
levels |
The home range level. |
tol |
Numeric scalar, indicating which difference of to stop. |
max.it |
Numeric scalar, indicating the maximum number of acceptable iterations. |
This implementation uses a bisection algorithm to the find the smallest value
value for the kernel bandwidth within range
that produces an home-range
isopleth at level
consisting of n
polygons. Note, no difference is
is made between the two dimensions.
list
with the calculated bandwidth, exit status and the number of iteration.
Kie, John G. "A rule-based ad hoc method for selecting a bandwidth in kernel home-range analyses." Animal Biotelemetry 1.1 (2013): 1-12.
Sometimes the percentage overlap between a spatial polygon an a home-range is required.
hr_overlap_feature(x, sf, direction = "hr_with_feature", feature_names = NULL)
hr_overlap_feature(x, sf, direction = "hr_with_feature", feature_names = NULL)
x |
A home-range estimate. |
sf |
An object of class |
direction |
The direction. |
feature_names |
optional feature names |
A tibble
Methods to calculate the overlap of two or more home-range estimates.
hr_overlap(x, ...) ## S3 method for class 'hr' hr_overlap(x, y, type = "hr", conditional = FALSE, ...) ## S3 method for class 'list' hr_overlap( x, type = "hr", conditional = FALSE, which = "consecutive", labels = NULL, ... )
hr_overlap(x, ...) ## S3 method for class 'hr' hr_overlap(x, y, type = "hr", conditional = FALSE, ...) ## S3 method for class 'list' hr_overlap( x, type = "hr", conditional = FALSE, which = "consecutive", labels = NULL, ... )
x , y
|
A home-range estimate |
... |
Further arguments, none implemented. |
type |
|
conditional |
|
which |
|
labels |
|
data.frame
with the isopleth level and area in units of the
coordinate reference system.
sfc
Convert a list column with many home-range estimates to a tibble
with a geometry column (as used by the sf
-package).
hr_to_sf(x, ...) ## S3 method for class 'tbl_df' hr_to_sf(x, col, ...)
hr_to_sf(x, ...) ## S3 method for class 'tbl_df' hr_to_sf(x, col, ...)
x |
A |
... |
Additional columns that should be transferred to the new |
col |
The column where the home |
A data.frame
with a simple feature column (from the sf
) package.
data("amt_fisher") hr <- amt_fisher |> nest(data = -id) |> mutate(hr = map(data, hr_mcp), n = map_int(data, nrow)) |> hr_to_sf(hr, id, n) hr <- amt_fisher |> nest(data = -id) |> mutate(hr = map(data, hr_kde), n = map_int(data, nrow)) |> hr_to_sf(hr, id, n)
data("amt_fisher") hr <- amt_fisher |> nest(data = -id) |> mutate(hr = map(data, hr_mcp), n = map_int(data, nrow)) |> hr_to_sf(hr, id, n) hr <- amt_fisher |> nest(data = -id) |> mutate(hr = map(data, hr_kde), n = map_int(data, nrow)) |> hr_to_sf(hr, id, n)
Obtain the utilization distribution of a probabilistic home-range estimate
hr_ud(x, ...)
hr_ud(x, ...)
x |
|
... |
Further arguments, none implemented. |
SpatRaster
Provides a very basic interface to leaflet
and lets the user inspect relocations on an interactive map.
inspect(x, ...) ## S3 method for class 'track_xy' inspect(x, popup = NULL, cluster = TRUE, ...)
inspect(x, ...) ## S3 method for class 'track_xy' inspect(x, popup = NULL, cluster = TRUE, ...)
x |
|
... |
Further arguments, none implemented. |
popup |
|
cluster |
|
An interactive leaflet
map.
Important, x
requires a valid coordinate reference system.
leaflet::leaflet()
data(sh) x <- track(x = sh$x, y = sh$y, crs = 31467) inspect(x) inspect(x, cluster = FALSE) inspect(x, popup = 1:nrow(x), cluster = FALSE)
data(sh) x <- track(x = sh$x, y = sh$y, crs = 31467) inspect(x) inspect(x, cluster = FALSE) inspect(x, popup = 1:nrow(x), cluster = FALSE)
Creates a formula without stratum variable
issf_drop_stratum(object, l)
issf_drop_stratum(object, l)
object |
|
l |
|
Creates a formula without movement variables
issf_w_form(object, l)
issf_w_form(object, l)
object |
|
l |
|
Calculate log-RSS(x1, x2) for a fitted RSF or (i)SSF
log_rss(object, ...) ## S3 method for class 'glm' log_rss(object, x1, x2, ci = NA, ci_level = 0.95, n_boot = 1000, ...) ## S3 method for class 'fit_clogit' log_rss(object, x1, x2, ci = NA, ci_level = 0.95, n_boot = 1000, ...)
log_rss(object, ...) ## S3 method for class 'glm' log_rss(object, x1, x2, ci = NA, ci_level = 0.95, n_boot = 1000, ...) ## S3 method for class 'fit_clogit' log_rss(object, x1, x2, ci = NA, ci_level = 0.95, n_boot = 1000, ...)
object |
|
... |
Further arguments, none implemented. |
x1 |
|
x2 |
|
ci |
|
ci_level |
|
n_boot |
|
This function assumes that the user would like to compare relative
selection strengths from at least one proposed location (x1
) to exactly
one reference location (x2
).
The objects object$model
, x1
, and x2
will be passed to
predict()
. Therefore, the columns of x1
and x2
must match
the terms in the model formula exactly.
Returns a list
of class log_rss
with four entries:
df
: A data.frame
with the covariates and the log_rss
x1
: A data.frame
with covariate values for x1
.
x2
: A data.frame
with covariate values for x2
.
formula
: The formula used to fit the model.
Brian J. Smith
Avgar, T., Lele, S.R., Keim, J.L., and Boyce, M.S.. (2017). Relative Selection Strength: Quantifying effect size in habitat- and step-selection inference. Ecology and Evolution, 7, 5322–5330.
Fieberg, J., Signer, J., Smith, B., & Avgar, T. (2021). A "How to" guide for interpreting parameters in habitat-selection analyses. Journal of Animal Ecology, 90(5), 1027-1043.
See Avgar et al. 2017 for details about relative selection strength.
Default plotting method available: plot.log_rss()
# RSF ------------------------------------------------------- # Fit an RSF, then calculate log-RSS to visualize results. # Load packages library(ggplot2) # Load data data("amt_fisher") amt_fisher_covar <- get_amt_fisher_covars() # Prepare data for RSF rsf_data <- amt_fisher |> filter(name == "Lupe") |> make_track(x_, y_, t_) |> random_points() |> extract_covariates(amt_fisher_covar$elevation) |> extract_covariates(amt_fisher_covar$popden) |> extract_covariates(amt_fisher_covar$landuse) |> mutate(lu = factor(landuse)) # Fit RSF m1 <- rsf_data |> fit_rsf(case_ ~ lu + elevation + popden) # Calculate log-RSS # data.frame of x1s x1 <- data.frame(lu = factor(50, levels = levels(rsf_data$lu)), elevation = seq(90, 120, length.out = 100), popden = mean(rsf_data$popden)) # data.frame of x2 (note factor levels should be same as model data) x2 <- data.frame(lu = factor(50, levels = levels(rsf_data$lu)), elevation = mean(rsf_data$elevation), popden = mean(rsf_data$popden)) # Calculate (use se method for confidence interval) logRSS <- log_rss(object = m1, x1 = x1, x2 = x2, ci = "se") # Plot ggplot(logRSS$df, aes(x = elevation_x1, y = log_rss)) + geom_hline(yintercept = 0, linetype = "dashed", color = "gray") + geom_ribbon(aes(ymin = lwr, ymax = upr), fill = "gray80") + geom_line() + xlab(expression("Elevation " * (x[1]))) + ylab("log-RSS") + ggtitle(expression("log-RSS" * (x[1] * ", " * x[2]))) + theme_bw() # SSF ------------------------------------------------------- # Fit an SSF, then calculate log-RSS to visualize results. # Load data data(deer) sh_forest <- get_sh_forest() # Prepare data for SSF ssf_data <- deer |> steps_by_burst() |> random_steps(n = 15) |> extract_covariates(sh_forest) |> mutate(forest = factor(forest, levels = 1:0, labels = c("forest", "non-forest")), cos_ta = cos(ta_), log_sl = log(sl_)) # Fit an SSF (note model = TRUE necessary for predict() to work) m2 <- ssf_data |> fit_clogit(case_ ~ forest + strata(step_id_), model = TRUE) # Calculate log-RSS # data.frame of x1s x1 <- data.frame(forest = factor(c("forest", "non-forest"))) # data.frame of x2 x2 <- data.frame(forest = factor("forest", levels = levels(ssf_data$forest))) # Calculate logRSS <- log_rss(object = m2, x1 = x1, x2 = x2, ci = "se") # Plot ggplot(logRSS$df, aes(x = forest_x1, y = log_rss)) + geom_hline(yintercept = 0, linetype = "dashed", color = "gray") + geom_errorbar(aes(ymin = lwr, ymax = upr), width = 0.25) + geom_point(size = 3) + xlab(expression("Forest Cover " * (x[1]))) + ylab("log-RSS") + ggtitle(expression("log-RSS" * (x[1] * ", " * x[2]))) + theme_bw()
# RSF ------------------------------------------------------- # Fit an RSF, then calculate log-RSS to visualize results. # Load packages library(ggplot2) # Load data data("amt_fisher") amt_fisher_covar <- get_amt_fisher_covars() # Prepare data for RSF rsf_data <- amt_fisher |> filter(name == "Lupe") |> make_track(x_, y_, t_) |> random_points() |> extract_covariates(amt_fisher_covar$elevation) |> extract_covariates(amt_fisher_covar$popden) |> extract_covariates(amt_fisher_covar$landuse) |> mutate(lu = factor(landuse)) # Fit RSF m1 <- rsf_data |> fit_rsf(case_ ~ lu + elevation + popden) # Calculate log-RSS # data.frame of x1s x1 <- data.frame(lu = factor(50, levels = levels(rsf_data$lu)), elevation = seq(90, 120, length.out = 100), popden = mean(rsf_data$popden)) # data.frame of x2 (note factor levels should be same as model data) x2 <- data.frame(lu = factor(50, levels = levels(rsf_data$lu)), elevation = mean(rsf_data$elevation), popden = mean(rsf_data$popden)) # Calculate (use se method for confidence interval) logRSS <- log_rss(object = m1, x1 = x1, x2 = x2, ci = "se") # Plot ggplot(logRSS$df, aes(x = elevation_x1, y = log_rss)) + geom_hline(yintercept = 0, linetype = "dashed", color = "gray") + geom_ribbon(aes(ymin = lwr, ymax = upr), fill = "gray80") + geom_line() + xlab(expression("Elevation " * (x[1]))) + ylab("log-RSS") + ggtitle(expression("log-RSS" * (x[1] * ", " * x[2]))) + theme_bw() # SSF ------------------------------------------------------- # Fit an SSF, then calculate log-RSS to visualize results. # Load data data(deer) sh_forest <- get_sh_forest() # Prepare data for SSF ssf_data <- deer |> steps_by_burst() |> random_steps(n = 15) |> extract_covariates(sh_forest) |> mutate(forest = factor(forest, levels = 1:0, labels = c("forest", "non-forest")), cos_ta = cos(ta_), log_sl = log(sl_)) # Fit an SSF (note model = TRUE necessary for predict() to work) m2 <- ssf_data |> fit_clogit(case_ ~ forest + strata(step_id_), model = TRUE) # Calculate log-RSS # data.frame of x1s x1 <- data.frame(forest = factor(c("forest", "non-forest"))) # data.frame of x2 x2 <- data.frame(forest = factor("forest", levels = levels(ssf_data$forest))) # Calculate logRSS <- log_rss(object = m2, x1 = x1, x2 = x2, ci = "se") # Plot ggplot(logRSS$df, aes(x = forest_x1, y = log_rss)) + geom_hline(yintercept = 0, linetype = "dashed", color = "gray") + geom_errorbar(aes(ymin = lwr, ymax = upr), width = 0.25) + geom_point(size = 3) + xlab(expression("Forest Cover " * (x[1]))) + ylab("log-RSS") + ggtitle(expression("log-RSS" * (x[1] * ", " * x[2]))) + theme_bw()
issf
-model object from scratchIn order to simulate from an issf
a
make_issf_model( coefs = c(sl_ = 0), sl = make_exp_distr(), ta = make_unif_distr() )
make_issf_model( coefs = c(sl_ = 0), sl = make_exp_distr(), ta = make_unif_distr() )
coefs |
A named vector with the coefficient values. |
sl |
The tentative step-length distribution. |
ta |
The tentative turn-angle distribution. |
An object of fit_clogit
.
An initial step for simulations. This step can either be created by defining a step from scratch or by using an observed step.
make_start(x, ...) ## S3 method for class 'numeric' make_start( x = c(0, 0), ta_ = 0, time = Sys.time(), dt = hours(1), crs = NA, ... ) ## S3 method for class 'track_xyt' make_start(x, ta_ = 0, dt = hours(1), ...) ## S3 method for class 'steps_xyt' make_start(x, ...)
make_start(x, ...) ## S3 method for class 'numeric' make_start( x = c(0, 0), ta_ = 0, time = Sys.time(), dt = hours(1), crs = NA, ... ) ## S3 method for class 'track_xyt' make_start(x, ta_ = 0, dt = hours(1), ...) ## S3 method for class 'steps_xyt' make_start(x, ...)
x |
|
... |
Further arguments, none implemented. |
ta_ |
|
time |
|
dt |
|
crs |
|
Functions to calculate metrics such as straightness, mean squared displacement (msd), intensity use,
sinuosity, mean turn angle correlation (tac
) of a track.
straightness(x, ...) cum_dist(x, ...) tot_dist(x, ...) msd(x, ...) intensity_use(x, ...) sinuosity(x, ...) tac(x, ...)
straightness(x, ...) cum_dist(x, ...) tot_dist(x, ...) msd(x, ...) intensity_use(x, ...) sinuosity(x, ...) tac(x, ...)
x |
|
... |
Further arguments, none implemented. |
The intensity use is calculated by dividing the total movement distance (tot_dist
) by the square of the area of movement (= minimum convex polygon 100).
A numeric vector of length one.
Abrahms B, Seidel DP, Dougherty E, Hazen EL, Bograd SJ, Wilson AM, McNutt JW, Costa DP, Blake S, Brashares JS, others (2017). “Suite of simple metrics reveals common movement syndromes across vertebrate taxa.” Movement ecology, 5(1), 12.
Almeida PJ, Vieira MV, Kajin M, Forero-Medina G, Cerqueira R (2010). “Indices of movement behaviour: conceptual background, effects of scale and location errors.” Zoologia (Curitiba), 27(5), 674–680.
Swihart RK, Slade NA (1985). “Testing for independence of observations in animal movements.” Ecology, 66(4), 1176–1184.
data(deer) tot_dist(deer) cum_dist(deer) straightness(deer) msd(deer) intensity_use(deer)
data(deer) tot_dist(deer) cum_dist(deer) straightness(deer) msd(deer) intensity_use(deer)
The function nsd()
calculates the net squared displacement (i.e., the squared distance from the first location of a track) for a track. The function add_nsd()
add a new column to a track or steps object with the nsd (the column name is nsd_
).
nsd(x, ...) ## S3 method for class 'track_xy' nsd(x, ...) add_nsd(x, ...) ## S3 method for class 'track_xy' add_nsd(x, ...) ## S3 method for class 'steps_xy' add_nsd(x, ...)
nsd(x, ...) ## S3 method for class 'track_xy' nsd(x, ...) add_nsd(x, ...) ## S3 method for class 'track_xy' add_nsd(x, ...) ## S3 method for class 'steps_xy' add_nsd(x, ...)
x |
|
... |
Further arguments, none implemented. |
Numeric vector (for nsd()
) and a tillbe with the original data with a new column for add_nsd()
.
od
is a wrapper around ctmm::occurrence
. See help(ctmm::occurrence)
for more details. rolling_od
estimates occurrence distributions for a subset of a track.
rolling_od(x, ...) ## S3 method for class 'track_xyt' rolling_od( x, trast, model = fit_ctmm(x, "bm"), res.space = 10, res.time = 10, n.points = 5, show.progress = TRUE, ... ) od(x, ...) ## S3 method for class 'track_xyt' od(x, trast, model = fit_ctmm(x, "bm"), res.space = 10, res.time = 10, ...)
rolling_od(x, ...) ## S3 method for class 'track_xyt' rolling_od( x, trast, model = fit_ctmm(x, "bm"), res.space = 10, res.time = 10, n.points = 5, show.progress = TRUE, ... ) od(x, ...) ## S3 method for class 'track_xyt' od(x, trast, model = fit_ctmm(x, "bm"), res.space = 10, res.time = 10, ...)
x |
|
... |
Further arguments, none implemented. |
trast |
|
model |
|
res.space |
|
res.time |
|
n.points |
|
show.progress |
|
Fleming, C. H., Fagan, W. F., Mueller, T., Olson, K. A., Leimgruber, P., & Calabrese, J. M. (2016). Estimating where and how animals travel: an optimal framework for path reconstruction from autocorrelated tracking data. Ecology.
data(deer) mini_deer <- deer[1:100, ] trast <- make_trast(mini_deer) md <- od(mini_deer, trast = trast) terra::plot(md) # rolling ud xx <- rolling_od(mini_deer, trast)
data(deer) mini_deer <- deer[1:100, ] trast <- make_trast(mini_deer) md <- od(mini_deer, trast = trast) terra::plot(md) # rolling ud xx <- rolling_od(mini_deer, trast)
Get parameters from a (fitted) distribution
sl_distr_params(x, ...) ## S3 method for class 'random_steps' sl_distr_params(x, ...) ## S3 method for class 'fit_clogit' sl_distr_params(x, ...) ta_distr_params(x, ...) ## S3 method for class 'random_steps' ta_distr_params(x, ...) ## S3 method for class 'fit_clogit' ta_distr_params(x, ...)
sl_distr_params(x, ...) ## S3 method for class 'random_steps' sl_distr_params(x, ...) ## S3 method for class 'fit_clogit' sl_distr_params(x, ...) ta_distr_params(x, ...) ## S3 method for class 'random_steps' ta_distr_params(x, ...) ## S3 method for class 'fit_clogit' ta_distr_params(x, ...)
x |
|
... |
None |
A list with the parameters of the distribution.
data(deer) d <- deer |> steps() |> random_steps() sl_distr_params(d) ta_distr_params(d)
data(deer) d <- deer |> steps() |> random_steps() sl_distr_params(d) ta_distr_params(d)
Plot step-length distribution
plot_sl(x, ...) ## S3 method for class 'fit_clogit' plot_sl(x, n = 1000, upper_quantile = 0.99, plot = TRUE, ...) ## S3 method for class 'random_steps' plot_sl(x, n = 1000, upper_quantile = 0.99, plot = TRUE, ...)
plot_sl(x, ...) ## S3 method for class 'fit_clogit' plot_sl(x, n = 1000, upper_quantile = 0.99, plot = TRUE, ...) ## S3 method for class 'random_steps' plot_sl(x, n = 1000, upper_quantile = 0.99, plot = TRUE, ...)
x |
|
... |
Further arguments, none implemented. |
n |
|
upper_quantile |
|
plot |
|
A plot of the step-length distribution.
data(deer) # with random steps deer[1:100, ] |> steps_by_burst() |> random_steps() |> plot_sl() deer[1:100, ] |> steps_by_burst() |> random_steps() |> plot_sl(upper_quantile = 0.5)
data(deer) # with random steps deer[1:100, ] |> steps_by_burst() |> random_steps() |> plot_sl() deer[1:100, ] |> steps_by_burst() |> random_steps() |> plot_sl(upper_quantile = 0.5)
Plot a home-range estimate
## S3 method for class 'hr' plot(x, add.relocations = TRUE, ...)
## S3 method for class 'hr' plot(x, add.relocations = TRUE, ...)
x |
A home-range estimate. |
add.relocations |
|
... |
Further arguments, none implemented. |
A plot
log_rss
objectDefault plot method for an object of class log_rss
## S3 method for class 'log_rss' plot(x, x_var1 = "guess", x_var2 = "guess", ...)
## S3 method for class 'log_rss' plot(x, x_var1 = "guess", x_var2 = "guess", ...)
x |
|
x_var1 |
|
x_var2 |
|
... |
|
This function provides defaults for a basic plot, but we encourage the user to carefully consider how to represent the patterns found in their habitat selection model.
The function log_rss()
is meant to accept a user-defined
input for x1
. The structure of x1
likely reflects how the user intended
to visualize the results. Therefore, it is possible to "guess" which covariate
the user would like to see on the x-axis by choosing the column from x1
with
the most unique values. Similarly, if there is a second column with multiple
unique values, that could be represented by a color. Note that if the user needs
to specify x_var1
, then we probably cannot guess x_var2
. Therefore, if the
user specifies x_var1 != "guess" & x_var2 == "guess"
, the function will return
an error.
This function uses integers to represent colors, and therefore the user can
change the default colors by specifying a custom palette()
before
calling the function.
A plot.
# Load data data("amt_fisher") amt_fisher_covar <- get_amt_fisher_covars() # Prepare data for RSF rsf_data <- amt_fisher |> filter(name == "Leroy") |> make_track(x_, y_, t_) |> random_points() |> extract_covariates(amt_fisher_covar$landuse) |> mutate(lu = factor(landuse)) # Fit RSF m1 <- rsf_data |> fit_rsf(case_ ~ lu) # Calculate log-RSS # data.frame of x1s x1 <- data.frame(lu = sort(unique(rsf_data$lu))) # data.frame of x2 (note factor levels should be same as model data) x2 <- data.frame(lu = factor(140, levels = levels(rsf_data$lu))) # Calculate logRSS <- log_rss(object = m1, x1 = x1, x2 = x2) # Plot plot(logRSS)
# Load data data("amt_fisher") amt_fisher_covar <- get_amt_fisher_covars() # Prepare data for RSF rsf_data <- amt_fisher |> filter(name == "Leroy") |> make_track(x_, y_, t_) |> random_points() |> extract_covariates(amt_fisher_covar$landuse) |> mutate(lu = factor(landuse)) # Fit RSF m1 <- rsf_data |> fit_rsf(case_ ~ lu) # Calculate log-RSS # data.frame of x1s x1 <- data.frame(lu = sort(unique(rsf_data$lu))) # data.frame of x2 (note factor levels should be same as model data) x2 <- data.frame(lu = factor(140, levels = levels(rsf_data$lu))) # Calculate logRSS <- log_rss(object = m1, x1 = x1, x2 = x2) # Plot plot(logRSS)
Plot an object of class uhc_data
## S3 method for class 'uhc_data' plot(x, ...)
## S3 method for class 'uhc_data' plot(x, ...)
x |
|
... |
Included for consistency with generic
|
Makes plots mimicking those in Fieberg et al. (2018), with the bootstrapped distribution in gray, the observed distribution in black, and the available distribution as a dashed red line.
Brian J. Smith
prep_uhc()
, conf_envelope()
,
plot.uhc_envelopes()
Plot an object of class uhc_envelopes
## S3 method for class 'uhc_envelopes' plot(x, ...)
## S3 method for class 'uhc_envelopes' plot(x, ...)
x |
|
... |
Included for consistency with generic
|
Makes plots mimicking those in Fieberg et al. (2018), with the
bootstrapped distribution in gray, the observed distribution in black,
and the available distribution as a dashed red line. This differs from
plot.uhc_data()
in that the bootstrapped distribution
(in gray) is drawn as a polygon rather than (many) lines, speeding up
plotting performance.
Brian J. Smith
prep_uhc()
, conf_envelope()
,
plot.uhc_data()
Creates data used to make used-habitat calibration plots
prep_uhc(object, test_dat, n_samp = 1000, n_dens = 512, verbose = TRUE) ## S3 method for class 'glm' prep_uhc(object, test_dat, n_samp = 1000, n_dens = 512, verbose = TRUE) ## S3 method for class 'fit_logit' prep_uhc(object, test_dat, n_samp = 1000, n_dens = 512, verbose = TRUE) ## S3 method for class 'fit_clogit' prep_uhc(object, test_dat, n_samp = 1000, n_dens = 512, verbose = TRUE)
prep_uhc(object, test_dat, n_samp = 1000, n_dens = 512, verbose = TRUE) ## S3 method for class 'glm' prep_uhc(object, test_dat, n_samp = 1000, n_dens = 512, verbose = TRUE) ## S3 method for class 'fit_logit' prep_uhc(object, test_dat, n_samp = 1000, n_dens = 512, verbose = TRUE) ## S3 method for class 'fit_clogit' prep_uhc(object, test_dat, n_samp = 1000, n_dens = 512, verbose = TRUE)
object |
|
test_dat |
|
n_samp |
|
n_dens |
|
verbose |
|
This function performs the heavy lifting of creating UHC plots.
It creates the data used later by the plot()
method, which actually
draws the UHC plots. This function (1) creates density plots of the used
and available locations from the test data, and (2) resamples the (a)
fitted coefficients and (b) test data (weighted by the exponential habitat
selection function) to create the distribution of used habitat under the
model.
Note that test_dat
should contain at least all of the variables that
appear in the model object
. Any further habitat variables in test_dat
will also have UHC plots generated, treating these variables as possible
candidate variables that are simply not included in this particular model.
Returns a list
of class uhc_data
with elements:
orig
: List of data.frame
s, one per variable (see vars
). Each
data.frame
contains the density plot data (x
and y
) for the original
used (dist == "U"
) and available (dist == "A"
) data.
samp
: List of data.frame
s, one per variable (see vars
). Each
data.frame
contains the density plot data (x
and y
) for each iteration
of bootstrap resampling (iter
).
vars
: Character vector with names of the habitat variables for which to
create UHC plots.
type
: Named character vector with the type for each of vars
(either
"numeric"
or "factor"
).
resp
: Character vector of length 1 with the name of the response
variable.
Brian J. Smith
Fieberg, J.R., Forester, J.D., Street, G.M., Johnson, D.H., ArchMiller, A.A., and Matthiopoulos, J. 2018. Used-habitat calibration plots: A new procedure for validating species distribution, resource selection, and step-selection models. Ecography 41:737–752.
See Fieberg et al. 2018 for details about UHC plots.
Default plotting method available: plot.uhc_data()
Coercion to data.frame
: as.data.frame.uhc_data()
Subsetting method: Extract.uhc_data
# Load packages library(amt) library(dplyr) library(terra) library(sf) # HSF ---------------------------------------------- # Load data data(uhc_hsf_locs) data(uhc_hab) hab <- rast(uhc_hab, type = "xyz", crs = "epsg:32612") # Convert "cover" layer to factor levels(hab[[4]]) <- data.frame(id = 1:3, cover = c("grass", "forest", "wetland")) # Split into train (80%) and test (20%) set.seed(1) uhc_hsf_locs$train <- rbinom(n = nrow(uhc_hsf_locs), size = 1, prob = 0.8) train <- uhc_hsf_locs[uhc_hsf_locs$train == 1, ] test <- uhc_hsf_locs[uhc_hsf_locs$train == 0, ] # Available locations avail_train <- random_points(st_as_sf(st_as_sfc(st_bbox(hab))), n = nrow(train) * 10) avail_test <- random_points(st_as_sf(st_as_sfc(st_bbox(hab))), n = nrow(test) * 10) # Combine with used train_dat <- train |> make_track(x, y, crs = 32612) |> mutate(case_ = TRUE) |> bind_rows(avail_train) |> # Attach covariates extract_covariates(hab) |> # Assign large weights to available mutate(weight = case_when( case_ ~ 1, !case_ ~ 5000 )) test_dat <- test |> make_track(x, y, crs = 32612) |> mutate(case_ = TRUE) |> bind_rows(avail_test) |> # Attach covariates extract_covariates(hab) |> # Assign large weights to available mutate(weight = case_when( case_ ~ 1, !case_ ~ 5000 )) # Fit (correct) HSF hsf1 <- glm(case_ ~ forage + temp + I(temp^2) + pred + cover, data = train_dat, family = binomial(), weights = weight) # Drop weights from 'test_dat' test_dat$weight <- NULL # Prep UHC plots uhc_dat <- prep_uhc(object = hsf1, test_dat = test_dat, n_samp = 500, verbose = TRUE) # Plot all variables plot(uhc_dat) # Plot only first variable plot(uhc_dat[1]) # Plot only "cover" variable plot(uhc_dat["cover"]) # Coerce to data.frame df <- as.data.frame(uhc_dat) # Simplify sampled lines to confidence envelopes conf <- conf_envelope(df) # Default plot for the envelopes version plot(conf)
# Load packages library(amt) library(dplyr) library(terra) library(sf) # HSF ---------------------------------------------- # Load data data(uhc_hsf_locs) data(uhc_hab) hab <- rast(uhc_hab, type = "xyz", crs = "epsg:32612") # Convert "cover" layer to factor levels(hab[[4]]) <- data.frame(id = 1:3, cover = c("grass", "forest", "wetland")) # Split into train (80%) and test (20%) set.seed(1) uhc_hsf_locs$train <- rbinom(n = nrow(uhc_hsf_locs), size = 1, prob = 0.8) train <- uhc_hsf_locs[uhc_hsf_locs$train == 1, ] test <- uhc_hsf_locs[uhc_hsf_locs$train == 0, ] # Available locations avail_train <- random_points(st_as_sf(st_as_sfc(st_bbox(hab))), n = nrow(train) * 10) avail_test <- random_points(st_as_sf(st_as_sfc(st_bbox(hab))), n = nrow(test) * 10) # Combine with used train_dat <- train |> make_track(x, y, crs = 32612) |> mutate(case_ = TRUE) |> bind_rows(avail_train) |> # Attach covariates extract_covariates(hab) |> # Assign large weights to available mutate(weight = case_when( case_ ~ 1, !case_ ~ 5000 )) test_dat <- test |> make_track(x, y, crs = 32612) |> mutate(case_ = TRUE) |> bind_rows(avail_test) |> # Attach covariates extract_covariates(hab) |> # Assign large weights to available mutate(weight = case_when( case_ ~ 1, !case_ ~ 5000 )) # Fit (correct) HSF hsf1 <- glm(case_ ~ forage + temp + I(temp^2) + pred + cover, data = train_dat, family = binomial(), weights = weight) # Drop weights from 'test_dat' test_dat$weight <- NULL # Prep UHC plots uhc_dat <- prep_uhc(object = hsf1, test_dat = test_dat, n_samp = 500, verbose = TRUE) # Plot all variables plot(uhc_dat) # Plot only first variable plot(uhc_dat[1]) # Plot only "cover" variable plot(uhc_dat["cover"]) # Coerce to data.frame df <- as.data.frame(uhc_dat) # Simplify sampled lines to confidence envelopes conf <- conf_envelope(df) # Default plot for the envelopes version plot(conf)
Sample random numbers from a distribution created within the amt
package.
random_numbers(x, n = 100, ...)
random_numbers(x, n = 100, ...)
x |
|
n |
|
... |
none implemented. |
A numeric vector.
Functions to generate random points within an animals home range. This is usually the first step for investigating habitat selection via Resource Selection Functions (RSF).
random_points(x, ...) ## S3 method for class 'hr' random_points(x, n = 100, type = "random", presence = NULL, ...) ## S3 method for class 'sf' random_points(x, n = 100, type = "random", presence = NULL, ...) ## S3 method for class 'track_xy' random_points(x, level = 1, hr = "mcp", n = nrow(x) * 10, type = "random", ...)
random_points(x, ...) ## S3 method for class 'hr' random_points(x, n = 100, type = "random", presence = NULL, ...) ## S3 method for class 'sf' random_points(x, n = 100, type = "random", presence = NULL, ...) ## S3 method for class 'track_xy' random_points(x, level = 1, hr = "mcp", n = nrow(x) * 10, type = "random", ...)
x |
|
... |
|
n |
|
type |
|
presence |
|
level |
|
hr |
|
A tibble
with the observed and random points and a new column case_
that indicates if a point is observed (case_ = TRUE
) or random (case_ TRUE
).
For objects of class track_xyt
the timestamp (t_
) is lost.
data(deer) # track_xyt --------------------------------------------------------------- # Default settings rp1 <- random_points(deer) plot(rp1) # Ten random points for each observed point rp <- random_points(deer, n = nrow(deer) * 10) plot(rp) # Within a home range ----------------------------------------------------- hr <- hr_mcp(deer, level = 1) # 100 random point within the home range rp <- random_points(hr, n = 100) plot(rp) # 100 regular point within the home range rp <- random_points(hr, n = 100, type = "regular") plot(rp) # 100 hexagonal point within the home range rp <- random_points(hr, n = 100, type = "hexagonal") plot(rp)
data(deer) # track_xyt --------------------------------------------------------------- # Default settings rp1 <- random_points(deer) plot(rp1) # Ten random points for each observed point rp <- random_points(deer, n = nrow(deer) * 10) plot(rp) # Within a home range ----------------------------------------------------- hr <- hr_mcp(deer, level = 1) # 100 random point within the home range rp <- random_points(hr, n = 100) plot(rp) # 100 regular point within the home range rp <- random_points(hr, n = 100, type = "regular") plot(rp) # 100 hexagonal point within the home range rp <- random_points(hr, n = 100, type = "hexagonal") plot(rp)
Function to generate a given number of random steps for each observed step.
random_steps(x, ...) ## S3 method for class 'numeric' random_steps( x, n_control = 10, angle = 0, rand_sl = random_numbers(make_exp_distr(), n = 1e+05), rand_ta = random_numbers(make_unif_distr(), n = 1e+05), ... ) ## S3 method for class 'steps_xy' random_steps( x, n_control = 10, sl_distr = fit_distr(x$sl_, "gamma"), ta_distr = fit_distr(x$ta_, "vonmises"), rand_sl = random_numbers(sl_distr, n = 1e+05), rand_ta = random_numbers(ta_distr, n = 1e+05), include_observed = TRUE, start_id = 1, ... ) ## S3 method for class 'bursted_steps_xyt' random_steps( x, n_control = 10, sl_distr = fit_distr(x$sl_, "gamma"), ta_distr = fit_distr(x$ta_, "vonmises"), rand_sl = random_numbers(sl_distr, n = 1e+05), rand_ta = random_numbers(ta_distr, n = 1e+05), include_observed = TRUE, ... )
random_steps(x, ...) ## S3 method for class 'numeric' random_steps( x, n_control = 10, angle = 0, rand_sl = random_numbers(make_exp_distr(), n = 1e+05), rand_ta = random_numbers(make_unif_distr(), n = 1e+05), ... ) ## S3 method for class 'steps_xy' random_steps( x, n_control = 10, sl_distr = fit_distr(x$sl_, "gamma"), ta_distr = fit_distr(x$ta_, "vonmises"), rand_sl = random_numbers(sl_distr, n = 1e+05), rand_ta = random_numbers(ta_distr, n = 1e+05), include_observed = TRUE, start_id = 1, ... ) ## S3 method for class 'bursted_steps_xyt' random_steps( x, n_control = 10, sl_distr = fit_distr(x$sl_, "gamma"), ta_distr = fit_distr(x$ta_, "vonmises"), rand_sl = random_numbers(sl_distr, n = 1e+05), rand_ta = random_numbers(ta_distr, n = 1e+05), include_observed = TRUE, ... )
x |
Steps. |
... |
Further arguments, none implemented. |
n_control |
|
angle |
|
rand_sl |
|
rand_ta |
|
sl_distr |
|
ta_distr |
|
include_observed |
|
start_id |
integer Index where the numbering for step ids start. |
A tibble
of class random_steps.
Simulate from an ssf model
random_steps_simple(start, sl_model, ta_model, n.control)
random_steps_simple(start, sl_model, ta_model, n.control)
start |
First step |
sl_model |
Step length model to use |
ta_model |
Turning angle model to use |
n.control |
How many alternative steps are considered each step |
Simulated trajectory
The range that in either x
, y
or both
directions, that a track covers.
range_x(x, ...) ## S3 method for class 'track_xy' range_x(x, ...) range_y(x, ...) ## S3 method for class 'track_xy' range_y(x, ...) range_both(x, ...) ## S3 method for class 'track_xy' range_both(x, ...)
range_x(x, ...) ## S3 method for class 'track_xy' range_x(x, ...) range_y(x, ...) ## S3 method for class 'track_xy' range_y(x, ...) range_both(x, ...) ## S3 method for class 'track_xy' range_both(x, ...)
x |
|
... |
Further arguments, none implemented. |
Numeric vector with the range.
From a fitted integrated step-selection function for a given position a redistribution kernel is calculated (i.e., the product of the movement kernel and the selection function).
redistribution_kernel( x = make_issf_model(), start = make_start(), map, fun = function(xy, map) { extract_covariates(xy, map, where = "both") }, covars = NULL, max.dist = get_max_dist(x), n.control = 1e+06, n.sample = 1, landscape = "continuous", compensate.movement = landscape == "discrete", normalize = TRUE, interpolate = FALSE, as.rast = FALSE, tolerance.outside = 0 )
redistribution_kernel( x = make_issf_model(), start = make_start(), map, fun = function(xy, map) { extract_covariates(xy, map, where = "both") }, covars = NULL, max.dist = get_max_dist(x), n.control = 1e+06, n.sample = 1, landscape = "continuous", compensate.movement = landscape == "discrete", normalize = TRUE, interpolate = FALSE, as.rast = FALSE, tolerance.outside = 0 )
x |
|
start |
|
map |
|
fun |
|
covars |
|
max.dist |
|
n.control |
|
n.sample |
|
landscape |
|
compensate.movement |
|
normalize |
|
interpolate |
|
as.rast |
|
tolerance.outside |
|
Removing relocations at the beginning and/or end of a track, that fall within a user specified period.
remove_capture_effect(x, ...) ## S3 method for class 'track_xyt' remove_capture_effect(x, start, end, ...)
remove_capture_effect(x, ...) ## S3 method for class 'track_xyt' remove_capture_effect(x, start, end, ...)
x |
An object of class |
... |
Further arguments, none implemented. |
start |
A |
end |
A |
A tibble
without observations that fall within the period of the capture effect.
Remove strata with missing values for observed steps
remove_incomplete_strata(x, ...) ## S3 method for class 'random_steps' remove_incomplete_strata(x, col = "ta_", ...)
remove_incomplete_strata(x, ...) ## S3 method for class 'random_steps' remove_incomplete_strata(x, col = "ta_", ...)
x |
An object of class |
... |
Further arguments, none implemented. |
col |
A character with the column name that will be scanned for missing values. |
An object of class random_steps
, where observed steps (case_ == TRUE
) with missing values (NA
) in the column col
are removed (including all random steps).
mini_deer <- deer[1:4, ] # The first step is removed, because we have `NA` turn angles. mini_deer |> steps() |> random_steps() |> remove_incomplete_strata() |> select(case_, ta_, step_id_)
mini_deer <- deer[1:4, ] # The first step is removed, because we have `NA` turn angles. mini_deer |> steps() |> random_steps() |> remove_incomplete_strata() |> select(case_, ta_, step_id_)
Extracts sampling period from a track_xyt
object
sampling_period(x, ...)
sampling_period(x, ...)
x |
|
... |
Further arguments, none implemented. |
Calculates SDR for an object of certain class
sdr(x, time_unit = "secs", append_na = TRUE, ...) ## S3 method for class 'track_xyt' sdr(x, time_unit = "secs", append_na = TRUE, ...)
sdr(x, time_unit = "secs", append_na = TRUE, ...) ## S3 method for class 'track_xyt' sdr(x, time_unit = "secs", append_na = TRUE, ...)
x |
|
time_unit |
|
append_na |
|
... |
Further arguments, none implemented. |
time_unit
defaults to seconds because calculate_sdr()
returns SDR in m^2/s. We assume the user is also working in a projected
CRS with units in meters, thus we expect SDR in m^2/s to be most relevant.
Brian J. Smith and Johannes Signer
calculate_sdr()
, get_displacement()
1500 GPS fixes of one red deer from northern Germany.
sh
sh
A data frame with 1500 rows and 4 variables:
the x-coordinate
the y-coordinate
the day of the relocation
the hour of the relocation
Verein für Wildtierforschung Dresden und Göttingen e.V.
Forest cover for the home range of one red deer from northern Germany.
sh_forest
sh_forest
A SpatRast
other
forest
JRC
A. Pekkarinen, L. Reithmaier, P. Strobl (2007): Pan-European Forest/Non-Forest mapping with Landsat ETM+ and CORINE Land Cover 2000 data.
Function to simulate a movement trajectory (path) from a redistribution kernel.
simulate_path(x, ...) ## Default S3 method: simulate_path(x, ...) ## S3 method for class 'redistribution_kernel' simulate_path(x, n.steps = 100, start = x$args$start, verbose = FALSE, ...)
simulate_path(x, ...) ## Default S3 method: simulate_path(x, ...) ## S3 method for class 'redistribution_kernel' simulate_path(x, n.steps = 100, start = x$args$start, verbose = FALSE, ...)
x |
|
... |
Further arguments, none implemented. |
n.steps |
|
start |
|
verbose |
|
Calculates two indices (mean squared displacement and linearity) to test for site fidelity. Significance testing is done by permuting step lengths and drawing turning angles from a uniform distribution.
site_fidelity(x, ...) ## S3 method for class 'steps_xy' site_fidelity(x, n = 100, alpha = 0.05, ...)
site_fidelity(x, ...) ## S3 method for class 'steps_xy' site_fidelity(x, n = 100, alpha = 0.05, ...)
x |
A track |
... |
None implemented |
n |
Numeric scalar. The number of simulated trajectories. |
alpha |
Numeric scalar. The alpha value used for the bootstrapping. |
A list of length 4. msd_dat
and li_dat
is the mean square distance and linearity for the real date. msd_sim
and 'li_sim“ are the mean square distances and linearities for the simulated trajectories.
Spencer, S. R., Cameron, G. N., & Swihart, R. K. (1990). Operationally defining home range: temporal dependence exhibited by hispid cotton rats. Ecology, 1817-1822.
# real data data(deer) ds <- deer |> steps_by_burst() site_fidelity(ds)
# real data data(deer) ds <- deer |> steps_by_burst() site_fidelity(ds)
Obtain the speed of a track.
speed(x, ...) ## S3 method for class 'track_xyt' speed(x, append_na = TRUE, ...)
speed(x, ...) ## S3 method for class 'track_xyt' speed(x, append_na = TRUE, ...)
x |
A |
... |
Further arguments, none implemented. |
append_na |
|
[numeric]
The speed in m/s
.
clogit
formula and returns a formula without the strata
and the
left-hand sideTakes a clogit
formula and returns a formula without the strata
and the
left-hand side
ssf_formula(formula)
ssf_formula(formula)
formula |
A formula object |
f1 <- case_ ~ x1 * x2 + strata(step_id_) ssf_formula(f1)
f1 <- case_ ~ x1 * x2 + strata(step_id_) ssf_formula(f1)
Given a fitted ssf, and new location the weights for each location is calculated
ssf_weights(xy, object, compensate.movement = FALSE)
ssf_weights(xy, object, compensate.movement = FALSE)
xy |
The new locations. |
object |
The the fitted (i)SSF. |
compensate.movement |
Whether or not for the transformation from polar to Cartesian coordinates is corrected. |
step_lengths
can be use to calculate step lengths of a track. direction_abs
and direction_rel
calculate the absolute and relative direction of steps. steps
converts a track_xy*
from a point representation to a step representation and automatically calculates step lengths and relative turning angles.
direction_abs(x, ...) ## S3 method for class 'track_xy' direction_abs( x, full_circle = FALSE, zero_dir = "E", clockwise = FALSE, append_last = TRUE, lonlat = FALSE, ... ) direction_rel(x, ...) ## S3 method for class 'track_xy' direction_rel(x, lonlat = FALSE, append_last = TRUE, zero_dir = "E", ...) step_lengths(x, ...) ## S3 method for class 'track_xy' step_lengths(x, lonlat = FALSE, append_last = TRUE, ...) steps_by_burst(x, ...) ## S3 method for class 'track_xyt' steps_by_burst(x, lonlat = FALSE, keep_cols = NULL, ...) steps(x, ...) ## S3 method for class 'track_xy' steps(x, lonlat = FALSE, keep_cols = NULL, ...) ## S3 method for class 'track_xyt' steps(x, lonlat = FALSE, keep_cols = NULL, diff_time_units = "auto", ...)
direction_abs(x, ...) ## S3 method for class 'track_xy' direction_abs( x, full_circle = FALSE, zero_dir = "E", clockwise = FALSE, append_last = TRUE, lonlat = FALSE, ... ) direction_rel(x, ...) ## S3 method for class 'track_xy' direction_rel(x, lonlat = FALSE, append_last = TRUE, zero_dir = "E", ...) step_lengths(x, ...) ## S3 method for class 'track_xy' step_lengths(x, lonlat = FALSE, append_last = TRUE, ...) steps_by_burst(x, ...) ## S3 method for class 'track_xyt' steps_by_burst(x, lonlat = FALSE, keep_cols = NULL, ...) steps(x, ...) ## S3 method for class 'track_xy' steps(x, lonlat = FALSE, keep_cols = NULL, ...) ## S3 method for class 'track_xyt' steps(x, lonlat = FALSE, keep_cols = NULL, diff_time_units = "auto", ...)
x |
|
... |
Further arguments, none implemented |
full_circle |
|
zero_dir |
|
clockwise |
|
append_last |
|
lonlat |
|
keep_cols |
|
diff_time_units |
|
dierctions_*()
returns NA
for 0 step lengths.
step_lengths
calculates the step lengths between points a long the path. The last value returned is NA
, because no observed step is 'started' at the last point. If lonlat = TRUE
, step_lengths()
wraps sf::st_distance()
.
[numeric]
For step_lengths()
and direction_*
a numeric vector. [data.frame]
For steps
and steps_by_burst
, containing the steps.
xy <- tibble( x = c(1, 4, 8, 8, 12, 12, 8, 0, 0, 4, 2), y = c(0, 0, 0, 8, 12, 12, 12, 12, 8, 4, 2)) trk <- make_track(xy, x, y) # append last direction_abs(trk, append_last = TRUE) direction_abs(trk, append_last = FALSE) # degrees direction_abs(trk) |> as_degree() # full circle or not: check direction_abs(trk, full_circle = TRUE) direction_abs(trk, full_circle = FALSE) direction_abs(trk, full_circle = TRUE) |> as_degree() direction_abs(trk, full_circle = FALSE) |> as_degree() # direction of 0 direction_abs(trk, full_circle = TRUE, zero_dir = "N") direction_abs(trk, full_circle = TRUE, zero_dir = "E") direction_abs(trk, full_circle = TRUE, zero_dir = "S") direction_abs(trk, full_circle = TRUE, zero_dir = "W") # clockwise or not direction_abs(trk, full_circle = TRUE, zero_dir = "N", clockwise = FALSE) direction_abs(trk, full_circle = TRUE, zero_dir = "N", clockwise = TRUE) # Bearing (i.e. azimuth): only for lon/lat direction_abs(trk, full_circle = FALSE, zero_dir = "N", lonlat = FALSE, clockwise = TRUE) direction_abs(trk, full_circle = FALSE, zero_dir = "N", lonlat = TRUE, clockwise = TRUE)
xy <- tibble( x = c(1, 4, 8, 8, 12, 12, 8, 0, 0, 4, 2), y = c(0, 0, 0, 8, 12, 12, 12, 12, 8, 4, 2)) trk <- make_track(xy, x, y) # append last direction_abs(trk, append_last = TRUE) direction_abs(trk, append_last = FALSE) # degrees direction_abs(trk) |> as_degree() # full circle or not: check direction_abs(trk, full_circle = TRUE) direction_abs(trk, full_circle = FALSE) direction_abs(trk, full_circle = TRUE) |> as_degree() direction_abs(trk, full_circle = FALSE) |> as_degree() # direction of 0 direction_abs(trk, full_circle = TRUE, zero_dir = "N") direction_abs(trk, full_circle = TRUE, zero_dir = "E") direction_abs(trk, full_circle = TRUE, zero_dir = "S") direction_abs(trk, full_circle = TRUE, zero_dir = "W") # clockwise or not direction_abs(trk, full_circle = TRUE, zero_dir = "N", clockwise = FALSE) direction_abs(trk, full_circle = TRUE, zero_dir = "N", clockwise = TRUE) # Bearing (i.e. azimuth): only for lon/lat direction_abs(trk, full_circle = FALSE, zero_dir = "N", lonlat = FALSE, clockwise = TRUE) direction_abs(trk, full_circle = FALSE, zero_dir = "N", lonlat = TRUE, clockwise = TRUE)
Returns a summary of sampling rates
summarize_sampling_rate(x, ...) ## S3 method for class 'track_xyt' summarize_sampling_rate( x, time_unit = "auto", summarize = TRUE, as_tibble = TRUE, ... ) summarize_sampling_rate_many(x, ...) ## S3 method for class 'track_xyt' summarize_sampling_rate_many(x, cols, time_unit = "auto", ...)
summarize_sampling_rate(x, ...) ## S3 method for class 'track_xyt' summarize_sampling_rate( x, time_unit = "auto", summarize = TRUE, as_tibble = TRUE, ... ) summarize_sampling_rate_many(x, ...) ## S3 method for class 'track_xyt' summarize_sampling_rate_many(x, cols, time_unit = "auto", ...)
x |
A |
... |
Further arguments, none implemented. |
time_unit |
|
summarize |
A logical. If |
as_tibble |
A logical. Should result be returned as |
cols |
Columns used for grouping. |
Depending on summarize
and as_tibble
, a vector, table or tibble.
data(deer) amt::summarize_sampling_rate(deer) data(amt_fisher) # Add the month amt_fisher |> mutate(yday = lubridate::yday(t_)) |> summarize_sampling_rate_many(c("id", "yday"))
data(deer) amt::summarize_sampling_rate(deer) data(amt_fisher) # Add the month amt_fisher |> mutate(yday = lubridate::yday(t_)) |> summarize_sampling_rate_many(c("id", "yday"))
Summarizes step lengths for a track_xy*
object
summarize_sl(x, ...)
summarize_sl(x, ...)
x |
|
... |
Further arguments, none implemented. |
Summarizes speeds for a track_xyt
object
summarize_speed(x, ...)
summarize_speed(x, ...)
x |
|
... |
Further arguments, none implemented. |
A convenience wrapper around suncalc::getSunlightTimes
to annotate if a fix was taken during day or night (optionally also include dawn and dusk).
time_of_day(x, ...) ## S3 method for class 'track_xyt' time_of_day(x, include.crepuscule = FALSE, ...) ## S3 method for class 'steps_xyt' time_of_day(x, include.crepuscule = FALSE, where = "end", ...)
time_of_day(x, ...) ## S3 method for class 'track_xyt' time_of_day(x, include.crepuscule = FALSE, ...) ## S3 method for class 'steps_xyt' time_of_day(x, include.crepuscule = FALSE, where = "end", ...)
x |
|
... |
Further arguments, none implemented. |
include.crepuscule |
|
where |
|
A tibble
with an additional column tod_
that contains the time of the day for each relocation.
data(deer) deer |> time_of_day() deer |> steps_by_burst() |> time_of_day() deer |> steps_by_burst() |> time_of_day(where = "start") deer |> steps_by_burst() |> time_of_day(where = "end") deer |> steps_by_burst() |> time_of_day(where = "both")
data(deer) deer |> time_of_day() deer |> steps_by_burst() |> time_of_day() deer |> steps_by_burst() |> time_of_day(where = "start") deer |> steps_by_burst() |> time_of_day(where = "end") deer |> steps_by_burst() |> time_of_day(where = "both")
track_*
Constructor to crate a track, the basic building block of the amt
package. A
track
is usually created from a set of x
and y
coordinates, possibly
time stamps, and any number of optional columns, such as id, sex, age, etc.
mk_track( tbl, .x, .y, .t, ..., crs = NA_crs_, order_by_ts = TRUE, check_duplicates = FALSE, all_cols = FALSE, verbose = FALSE ) make_track( tbl, .x, .y, .t, ..., crs = NA_crs_, order_by_ts = TRUE, check_duplicates = FALSE, all_cols = FALSE, verbose = FALSE ) track(x, y, t, ..., crs = NULL)
mk_track( tbl, .x, .y, .t, ..., crs = NA_crs_, order_by_ts = TRUE, check_duplicates = FALSE, all_cols = FALSE, verbose = FALSE ) make_track( tbl, .x, .y, .t, ..., crs = NA_crs_, order_by_ts = TRUE, check_duplicates = FALSE, all_cols = FALSE, verbose = FALSE ) track(x, y, t, ..., crs = NULL)
tbl |
|
.x , .y , .t
|
|
... |
|
crs |
|
order_by_ts |
|
check_duplicates |
|
all_cols |
|
verbose |
|
x , y
|
|
t |
|
If t
was provided an object of class track_xyt
is returned
otherwise a track_xy
.
Functions to only selects relocations that can be aligned with a new time series (within some tolerance).
track_align(x, ...) ## S3 method for class 'track_xyt' track_align(x, new.times, tolerance, ...)
track_align(x, ...) ## S3 method for class 'track_xyt' track_align(x, new.times, tolerance, ...)
x |
A track. |
... |
Further arguments, none implemented. |
new.times |
The new time trajectory. |
tolerance |
The tolerance between observed time stamps and new time stamps in seconds. This should be either a vector of length 1 or length |
A track_xyt
.
Function to resample a track at a predefined sampling rate within some tolerance.
track_resample(x, ...) ## S3 method for class 'track_xyt' track_resample(x, rate = hours(2), tolerance = minutes(15), start = 1, ...)
track_resample(x, ...) ## S3 method for class 'track_xyt' track_resample(x, rate = hours(2), tolerance = minutes(15), start = 1, ...)
x |
A |
... |
Further arguments, none implemented. |
rate |
A lubridate |
tolerance |
A lubridate |
start |
A integer scalar, that gives the relocation at which the sampling rate starts. |
A resampled track_xyt
.
Subsets a track_xyt
object
tracked_from_to(x, from = min(x$t_), to = max(x$t_))
tracked_from_to(x, from = min(x$t_), to = max(x$t_))
x |
|
from |
|
to |
|
Transforms the CRS for a track.
transform_coords(x, ...) ## S3 method for class 'track_xy' transform_coords(x, crs_to, crs_from, ...) transform_crs(x, ...)
transform_coords(x, ...) ## S3 method for class 'track_xy' transform_coords(x, crs_to, crs_from, ...) transform_crs(x, ...)
x |
|
... |
Further arguments, none implemented. |
crs_to |
|
crs_from |
|
A track with transformed coordinates.
sf::st_transform
data(deer) get_crs(deer) # project to geographical coordinates (note the CRS is taken automatically from the object deer). d1 <- transform_coords(deer, crs_to = 4326)
data(deer) get_crs(deer) # project to geographical coordinates (note the CRS is taken automatically from the object deer). d1 <- transform_coords(deer, crs_to = 4326)
For some home-range estimation methods (e.g., KDE) a template raster is needed. This functions helps to quickly create such a template raster.
make_trast(x, ...) ## S3 method for class 'track_xy' make_trast(x, factor = 1.5, res = max(c(extent_max(x)/100, 1e-09)), ...)
make_trast(x, ...) ## S3 method for class 'track_xy' make_trast(x, factor = 1.5, res = max(c(extent_max(x)/100, 1e-09)), ...)
x |
|
... |
Further arguments, none implemented. |
factor |
|
res |
|
A RastLayer
without values.
Internal function to summarize distribution of numeric or factor variables
ua_distr(name, type, data, lims, resp, n_dens, avail = TRUE)
ua_distr(name, type, data, lims, resp, n_dens, avail = TRUE)
name |
|
type |
|
data |
|
lims |
|
resp |
|
n_dens |
|
avail |
|
Simulated habitat rasters for demonstrating UHC plots
uhc_hab
uhc_hab
A RasterStack
with 1600 cells and 7 variables:
Forage biomass in g/m^2^ (resource)
mean annual temperature in °C (condition)
predator density in predators/100 km^2^ (risk)
landcover (forest > grassland > wetland)
distance to the wetland landcover (no effect)
distance to the centroid of the raster (no effect)
random integers (no effect)
Simulated HSF location data for demonstrating UHC plots
uhc_hsf_locs
uhc_hsf_locs
A data.frame
with 2000 rows and 2 variables:
x-coordinate in UTM Zone 12 (EPSG: 32612)
Y-coordinate in UTM Zone 12 (EPSG: 32612)
These data were simulated assuming an ordinary habitat selection function (HSF), i.e., all points are independent rather than arising from an underlying movement model.
True parameter values are:
forage
= log(5)/500 (resource)
temp^2
= -1 * log(2)/36 (condition; quadratic term)
temp
= (log(2)/36) * 26 (condition; linear term)
pred
= log(0.25)/5 (risk)
cover == "forest"
= log(2) (grassland is intercept)
cover == "wetland"
= log(1/2) (grassland is intercept)
Note: temp
is modeled as a quadratic term, with the strongest selection
occurring at 13 °C and all other temperatures less selected.
Note: dist_to_water
, dist_to_cent
, and rand
have no real effect
on our animal's selection and are included for demonstration purposes.
Simulated iSSF location data for demonstrating UHC plots
uhc_issf_locs
uhc_issf_locs
A data.frame
with 371 rows and 3 variables:
x-coordinate in UTM Zone 12 (EPSG: 32612)
Y-coordinate in UTM Zone 12 (EPSG: 32612)
timestamp of location (timezone "US/Mountain")
These data were simulated assuming an movement model, i.e., iSSA.
True movement-free habitat selection parameter values are:
forage
= log(8)/500 (resource)
temp^2
= -1 * log(8)/36 (condition; quadratic term)
temp
= (log(8)/36) * 26 (condition; linear term)
pred
= log(0.2)/5 (risk)
cover == "forest"
= log(2) (grassland is intercept)
cover == "wetland"
= log(1/2) (grassland is intercept)
dist_to_cent
= -1 * log(10)/500 (keeps trajectory away from boundary)
Note: temp
is modeled as a quadratic term, with the strongest selection
occurring at 13 °C and all other temperatures less selected.
Note: dist_to_water
and rand
have no real effect
on our animal's selection and are included for demonstration purposes.
True selection-free movement distributions are:
Step length: gamma(shape = 3, scale = 25)
Turn angle: vonMises(mu = 0, kappa = 0.5)
amt_distr
Functions to update amt_distr
from iSSF coefficients
update_gamma(dist, beta_sl, beta_log_sl) update_exp(dist, beta_sl) update_hnorm(dist, beta_sl_sq) update_lnorm(dist, beta_log_sl, beta_log_sl_sq) update_vonmises(dist, beta_cos_ta)
update_gamma(dist, beta_sl, beta_log_sl) update_exp(dist, beta_sl) update_hnorm(dist, beta_sl_sq) update_lnorm(dist, beta_log_sl, beta_log_sl_sq) update_vonmises(dist, beta_cos_ta)
dist |
|
beta_sl |
|
beta_log_sl |
|
beta_sl_sq |
|
beta_log_sl_sq |
|
beta_cos_ta |
|
These functions are called internally by
update_sl_distr()
and update_ta_distr()
.
However, those simple functions assume that the selection-free step-length
and turn-angle distributions are constant (i.e., they do not depend on
covariates). In the case of interactions between movement parameters and
covariates, the user will want to manually access these functions to update
their selection-free movement distributions.
A distribution
# sh_forest <- get_sh_forest() # # Fit an SSF, then update movement parameters. # # #Prepare data for SSF # ssf_data <- deer |> # steps_by_burst() |> # random_steps(n = 15) |> # extract_covariates(sh_forest) |> # mutate(forest = factor(forest, levels = 1:0, # labels = c("forest", "non-forest")), # cos_ta_ = cos(ta_), # log_sl_ = log(sl_)) # # # Check tentative distributions # # Step length # attr(ssf_data, "sl_") # # Turning angle # attr(ssf_data, "ta_") # # # Fit an iSSF (note model = TRUE necessary for predict() to work) # m1 <- ssf_data |> # fit_issf(case_ ~ forest * (sl_ + log_sl_ + cos_ta_) + # strata(step_id_), model = TRUE) # # # Update forest step lengths (the reference level) # forest_sl <- update_gamma(m1$sl_, # beta_sl = m1$model$coefficients["sl_"], # beta_log_sl = m1$model$coefficients["log_sl_"]) # # # Update non-forest step lengths # nonforest_sl <- update_gamma(m1$sl_, # beta_sl = m1$model$coefficients["sl_"] + # m1$model$coefficients["forestnon-forest:sl_"], # beta_log_sl = m1$model$coefficients["log_sl_"] + # m1$model$coefficients["forestnon-forest:log_sl_"]) # # # Update forest turn angles (the reference level) # forest_ta <- update_vonmises(m1$ta_, # beta_cos_ta = m1$model$coefficients["cos_ta_"]) # # # Update non-forest turn angles # nonforest_ta <- update_vonmises(m1$ta_, # beta_cos_ta = m1$model$coefficients["cos_ta_"] + # m1$model$coefficients["forestnon-forest:cos_ta_"]) #
# sh_forest <- get_sh_forest() # # Fit an SSF, then update movement parameters. # # #Prepare data for SSF # ssf_data <- deer |> # steps_by_burst() |> # random_steps(n = 15) |> # extract_covariates(sh_forest) |> # mutate(forest = factor(forest, levels = 1:0, # labels = c("forest", "non-forest")), # cos_ta_ = cos(ta_), # log_sl_ = log(sl_)) # # # Check tentative distributions # # Step length # attr(ssf_data, "sl_") # # Turning angle # attr(ssf_data, "ta_") # # # Fit an iSSF (note model = TRUE necessary for predict() to work) # m1 <- ssf_data |> # fit_issf(case_ ~ forest * (sl_ + log_sl_ + cos_ta_) + # strata(step_id_), model = TRUE) # # # Update forest step lengths (the reference level) # forest_sl <- update_gamma(m1$sl_, # beta_sl = m1$model$coefficients["sl_"], # beta_log_sl = m1$model$coefficients["log_sl_"]) # # # Update non-forest step lengths # nonforest_sl <- update_gamma(m1$sl_, # beta_sl = m1$model$coefficients["sl_"] + # m1$model$coefficients["forestnon-forest:sl_"], # beta_log_sl = m1$model$coefficients["log_sl_"] + # m1$model$coefficients["forestnon-forest:log_sl_"]) # # # Update forest turn angles (the reference level) # forest_ta <- update_vonmises(m1$ta_, # beta_cos_ta = m1$model$coefficients["cos_ta_"]) # # # Update non-forest turn angles # nonforest_ta <- update_vonmises(m1$ta_, # beta_cos_ta = m1$model$coefficients["cos_ta_"] + # m1$model$coefficients["forestnon-forest:cos_ta_"]) #
Update tentative step length or turning angle distribution from a fitted iSSF.
update_sl_distr( object, beta_sl = "sl_", beta_log_sl = "log_sl_", beta_sl_sq = "sl_sq_", beta_log_sl_sq = "log_sl_sq_", ... ) update_ta_distr(object, beta_cos_ta = "cos_ta_", ...)
update_sl_distr( object, beta_sl = "sl_", beta_log_sl = "log_sl_", beta_sl_sq = "sl_sq_", beta_log_sl_sq = "log_sl_sq_", ... ) update_ta_distr(object, beta_cos_ta = "cos_ta_", ...)
object |
|
beta_sl |
|
beta_log_sl |
|
beta_sl_sq |
|
beta_log_sl_sq |
|
... |
Further arguments, none implemented. |
beta_cos_ta |
|
An amt_distr
object, which consists of a list with the name
of
the distribution and its parameters (saved in params
).
Brian J. Smith and Johannes Signer
Fieberg J, Signer J, Smith BJ, Avgar T (2020). “A “How-to” Guide for Interpreting Parameters in Resource-and Step-Selection Analyses.” bioRxiv.
Wrapper to fit a distribution to data fit_distr()
# Fit an SSF, then update movement parameters. data(deer) mini_deer <- deer[1:100, ] sh_forest <- get_sh_forest() # Prepare data for SSF ssf_data <- mini_deer |> steps_by_burst() |> random_steps(n = 15) |> extract_covariates(sh_forest) |> mutate(forest = factor(forest, levels = 1:0, labels = c("forest", "non-forest")), cos_ta_ = cos(ta_), log_sl_ = log(sl_)) # Check tentative distributions # Step length sl_distr_params(ssf_data) attr(ssf_data, "sl_") # Turning angle ta_distr_params(ssf_data) # Fit an iSSF m1 <- ssf_data |> fit_issf(case_ ~ forest + sl_ + log_sl_ + cos_ta_ + strata(step_id_)) # Update step length distribution new_gamma <- update_sl_distr(m1) # Update turning angle distribution new_vm <- update_ta_distr(m1) # It is also possible to use different step length distributions # exponential step-length distribution s2 <- mini_deer |> steps_by_burst() s2 <- random_steps(s2, sl_distr = fit_distr(s2$sl_, "exp")) m2 <- s2 |> fit_clogit(case_ ~ sl_ + strata(step_id_)) update_sl_distr(m2) # half normal step-length distribution s3 <- mini_deer |> steps_by_burst() s3 <- random_steps(s3, sl_distr = fit_distr(s3$sl_, "hnorm")) m3 <- s3 |> mutate(sl_sq_ = sl_^2) |> fit_clogit(case_ ~ sl_sq_ + strata(step_id_)) update_sl_distr(m3) # log normal step-length distribution s4 <- mini_deer |> steps_by_burst() s4 <- random_steps(s4, sl_distr = fit_distr(s4$sl_, "lnorm")) m4 <- s4 |> mutate(log_sl_ = log(sl_), log_sl_sq_ = log(sl_)^2) |> fit_clogit(case_ ~ log_sl_ + log_sl_sq_ + strata(step_id_)) update_sl_distr(m4)
# Fit an SSF, then update movement parameters. data(deer) mini_deer <- deer[1:100, ] sh_forest <- get_sh_forest() # Prepare data for SSF ssf_data <- mini_deer |> steps_by_burst() |> random_steps(n = 15) |> extract_covariates(sh_forest) |> mutate(forest = factor(forest, levels = 1:0, labels = c("forest", "non-forest")), cos_ta_ = cos(ta_), log_sl_ = log(sl_)) # Check tentative distributions # Step length sl_distr_params(ssf_data) attr(ssf_data, "sl_") # Turning angle ta_distr_params(ssf_data) # Fit an iSSF m1 <- ssf_data |> fit_issf(case_ ~ forest + sl_ + log_sl_ + cos_ta_ + strata(step_id_)) # Update step length distribution new_gamma <- update_sl_distr(m1) # Update turning angle distribution new_vm <- update_ta_distr(m1) # It is also possible to use different step length distributions # exponential step-length distribution s2 <- mini_deer |> steps_by_burst() s2 <- random_steps(s2, sl_distr = fit_distr(s2$sl_, "exp")) m2 <- s2 |> fit_clogit(case_ ~ sl_ + strata(step_id_)) update_sl_distr(m2) # half normal step-length distribution s3 <- mini_deer |> steps_by_burst() s3 <- random_steps(s3, sl_distr = fit_distr(s3$sl_, "hnorm")) m3 <- s3 |> mutate(sl_sq_ = sl_^2) |> fit_clogit(case_ ~ sl_sq_ + strata(step_id_)) update_sl_distr(m3) # log normal step-length distribution s4 <- mini_deer |> steps_by_burst() s4 <- random_steps(s4, sl_distr = fit_distr(s4$sl_, "lnorm")) m4 <- s4 |> mutate(log_sl_ = log(sl_), log_sl_sq_ = log(sl_)^2) |> fit_clogit(case_ ~ log_sl_ + log_sl_sq_ + strata(step_id_)) update_sl_distr(m4)