.acs | R Documentation |
This function is the back-end of the acoustic-container (AC) and acoustic-container depth-contour (ACDC) algorithms.
.acs(
acoustics,
archival = NULL,
step = 120,
round_ts = TRUE,
plot_ts = TRUE,
bathy,
map = NULL,
detection_containers,
detection_kernels = NULL,
detection_kernels_overlap = NULL,
detection_time_window = 5,
mobility,
calc_depth_error = function(...) matrix(c(-2.5, 2.5), nrow = 2),
normalise = TRUE,
chunk = 1L,
save_record_spatial = 1L,
write_record_spatial_for_pf = NULL,
save_args = TRUE,
verbose = TRUE,
con = "",
progress = 1L,
check = TRUE
)
acoustics |
A dataframe that contains passive acoustic telemetry detection time series for a specific individual (see |
archival |
For the ACDC algorithm, |
step |
A number that defines the time step length (s) between consecutive detections. If |
round_ts |
A logical input that defines whether or not to round the |
plot_ts |
A logical input that defines whether or not to the plot movement series before the algorithm is initiated. |
bathy |
A |
map |
(optional) A blank |
detection_containers |
A list of detection containers, with one element for each number from |
detection_kernels |
A named list of detection probability kernels, from |
detection_kernels_overlap |
(optional) A named list (the ‘list_by_receiver’ element from |
detection_time_window |
(optional) A number that defines the maximum duration (s) between consecutive detections at different receivers such that they can be said to have occurred at ‘effectively the same time’. This indicates that the same transmission was detected by multiple receivers. If |
mobility |
A number that defines the distance (m) that an individual could move in the time steps between acoustic detections (see also |
calc_depth_error |
In the ACDC algorithm, |
normalise |
A logical variable that defines whether or not to normalise the map of possible locations at each time step so that they sum to one. |
chunk |
An integer that defines the chunk ID (from |
save_record_spatial |
An integer vector that defines the time steps for which to save a record of the spatial information from each time step. |
write_record_spatial_for_pf |
(optional) A named list, passed to |
save_args |
A logical input that defines whether or not to include a list of function arguments in the outputs. This can be switched off if the function is applied iteratively. |
verbose |
A logical variable that defines whether or not to print messages to the console or to file to relay function progress. If |
con |
If |
progress |
(optional) If the algorithm is implemented step-wise, |
check |
A logical input that defines whether or not to check function inputs. This can be switched off to improve computation time when the function is applied iteratively or via a front-end function (e.g., |
The function returns an acdc_record-class
object with the following elements: ‘map’, ‘record’, ‘time’, ‘args’, ‘chunks’ and ‘simplify’. The main output element is the ‘map’ RasterLayer that shows where the individual could have spent more or less time over the duration of the movement time series. The ‘record’ element records time-specific information on the possible locations of the individual, and can be used to plot maps of specific time points or to produce animations (for the time steps specified by save_record_spatial
). The ‘time’ element is a dataframe that defines the times of sequential stages in the algorithm's progression, providing a record of computation time. The ‘args’ element is a named list of user inputs that record the parameters used to generate the outputs (if save_args = TRUE
, otherwise the ‘args’ element is NULL
).
Edward Lavender
The front-end functions ac
and acdc
call .acs_pl
which in turn calls this function. acs_setup_containers
defines the detection containers required by this function. acs_setup_mobility
is used to examine the assumption of the constant ‘mobility’ parameter. acs_setup_detection_kernels
produces detection probability kernels for incorporation into the function. For calls via ac
and acdc
, acdc_simplify
simplifies the outputs and acdc_plot_trace
, acdc_plot_record
and acdc_animate_record
visualise the results.
#### Step (1) Prepare study site grid
# Grid resolution needs to be sufficiently high to capture detection probability/movement
# And sufficiently low to minimise computational costs
blank <- raster::raster(raster::extent(dat_gebco), res = c(75, 75))
gebco <- raster::resample(dat_gebco, blank)
#### Step (2) Implement setup_acdc_*() steps
# ... Define detection containers required for algorithm(s) (see acs_setup_containers())
#### Step (3) Prepare movement time series for algorithm(s)
# Add required columns to dataframes:
dat_acoustics$timestamp_num <- as.numeric(dat_acoustics$timestamp)
dat_archival$timestamp_num <- as.numeric(dat_archival$timestamp)
# Focus on an example individual
id <- 25
acc <- dat_acoustics[dat_acoustics$individual_id == id, ]
arc <- dat_archival[dat_archival$individual_id == id, ]
# Focus on the subset of data for which we have both acoustic and archival detections
acc <- acc[acc$timestamp >= min(arc$timestamp) - 2 * 60 &
acc$timestamp <= max(arc$timestamp) + 2 * 60, ]
arc <- arc[arc$timestamp >= min(acc$timestamp) - 2 * 60 &
arc$timestamp <= max(acc$timestamp) + 2 * 60, ]
# We'll focus on a one day period with overlapping detection/depth time series for speed
end <- as.POSIXct("2016-03-18")
acc <- acc[acc$timestamp <= end, ]
arc <- arc[arc$timestamp <= end, ]
arc <- arc[arc$timestamp >= min(acc$timestamp) - 2 * 60 &
arc$timestamp <= max(acc$timestamp) + 2 * 60, ]
# Process time series
# ... Detections should be rounded to the nearest step
# ... Duplicate detections
# ... ... (of the same individual at the same receiver in the same step)
# ... ... should be dropped
# ... This step has already been implemented for dat_acoustics.
#### Example (1) Implement AC algorithm with default arguments
out_acdc_1 <- flapper:::.acs(
acoustics = acc,
bathy = gebco,
detection_containers = dat_containers,
mobility = 200
)
#### Example (2) Implement ACDC algorithm with default arguments
out_acdc_2 <- flapper:::.acs(
acoustics = acc,
archival = arc,
bathy = gebco,
detection_containers = dat_containers,
mobility = 200,
calc_depth_error = function(...) matrix(c(-2.5, 2.5), nrow = 2)
)
#### Example (3) Implement AC or ACDC algorithm with detection probability kernels
## (A) Get detection container overlaps
# Define receiver locations as a SpatialPointsDataFrame object with a UTM CRS
proj_wgs84 <- sp::CRS(SRS_string = "EPSG:4326")
proj_utm <- sp::CRS(SRS_string = "EPSG:32629")
rownames(dat_moorings) <- dat_moorings$receiver_id
xy <- sp::SpatialPoints(
dat_moorings[, c("receiver_long", "receiver_lat")],
proj_wgs84
)
xy <- sp::spTransform(xy, proj_utm)
xy <- sp::SpatialPointsDataFrame(xy, dat_moorings[, c(
"receiver_id",
"receiver_start_date",
"receiver_end_date"
)])
# Get detection overlap(s) as a SpatialPolygonsDataFrame
containers <- get_detection_containers(
xy = sp::SpatialPoints(xy),
detection_range = 425,
coastline = dat_coast,
byid = TRUE
)
containers_df <- dat_moorings[, c(
"receiver_id",
"receiver_start_date",
"receiver_end_date"
)]
row.names(containers_df) <- names(containers)
containers <- sp::SpatialPolygonsDataFrame(containers, containers_df)
overlaps <- get_detection_containers_overlap(containers = containers)
## (B) Define detection probability function based on distance and detection_range
calc_dpr <-
function(x) {
ifelse(x <= 425, stats::plogis(2.5 + -0.02 * x), 0)
}
## (C) Get detection kernels (a slow step)
kernels <- acs_setup_detection_kernels(
xy = xy,
containers = dat_containers,
overlaps = overlaps,
calc_detection_pr = calc_dpr,
bathy = gebco
)
## (D) Implement algorithm
out_acdc_3 <- flapper:::.acs(
acoustics = acc,
archival = arc,
bathy = gebco,
detection_containers = dat_containers,
detection_kernels = kernels,
detection_kernels_overlap = overlaps$list_by_receiver,
detection_time_window = 10,
mobility = 200,
calc_depth_error = function(...) matrix(c(-2.5, 2.5), nrow = 2)
)
#### Example (4): Compare outputs with/without detection probability and normalisation
## Without detection kernels
out_acdc_4a <- flapper:::.acs(
acoustics = acc,
archival = arc,
bathy = gebco,
detection_containers = dat_containers,
mobility = 200,
calc_depth_error = function(...) matrix(c(-2.5, 2.5), nrow = 2)
)
out_acdc_4b <- flapper:::.acs(
acoustics = acc,
archival = arc,
bathy = gebco,
detection_containers = dat_containers,
normalise = TRUE,
mobility = 200,
calc_depth_error = function(...) matrix(c(-2.5, 2.5), nrow = 2)
)
## With detection kernels
out_acdc_4c <- flapper:::.acs(
acoustics = acc,
archival = arc,
bathy = gebco,
detection_containers = dat_containers,
detection_kernels = kernels,
mobility = 200,
calc_depth_error = function(...) matrix(c(-2.5, 2.5), nrow = 2)
)
out_acdc_4d <- flapper:::.acs(
acoustics = acc,
archival = arc,
bathy = gebco,
detection_containers = dat_containers,
detection_kernels = kernels, normalise = TRUE,
mobility = 200,
calc_depth_error = function(...) matrix(c(-2.5, 2.5), nrow = 2)
)
## Comparison of final maps
pp <- par(mfrow = c(2, 2))
raster::plot(out_acdc_4a$map, main = "4a")
raster::plot(out_acdc_4b$map, main = "4b")
raster::plot(out_acdc_4c$map, main = "4c")
raster::plot(out_acdc_4d$map, main = "4d")
par(pp)
#### Example (5): Implement AC or ACDC algorithm and write messages to file via 'con'
out_acdc_5 <- flapper:::.acs(
acoustics = acc,
archival = arc,
bathy = gebco,
detection_containers = dat_containers,
mobility = 200,
calc_depth_error = function(...) matrix(c(-2.5, 2.5), nrow = 2),
verbose = TRUE,
con = paste0(tempdir(), "/", "acdc_log.txt")
)
# Check log
utils::head(readLines(paste0(tempdir(), "/", "acdc_log.txt")))
utils::tail(readLines(paste0(tempdir(), "/", "acdc_log.txt")))
#### Example (6): Implement AC or ACDC algorithm and return spatial information
# Specify save_record_spatial = NULL to include spatial information for all time steps
# ... (used for plotting) or a vector to include this information for specific time steps
out_acdc_6 <- flapper:::.acs(
acoustics = acc,
archival = arc,
bathy = gebco,
detection_containers = dat_containers,
mobility = 200,
calc_depth_error = function(...) matrix(c(-2.5, 2.5), nrow = 2),
save_record_spatial = NULL,
verbose = TRUE,
con = paste0(tempdir(), "/", "acdc_log.txt")
)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.