pf_kud_ | R Documentation |
These functions are wrappers designed to apply kernel utilisation distribution (KUD) estimation to the outputs of a particle filtering (PF) algorithm. To implement these routines, an (a) pf_archive-class
object from pf
(plus pf_simplify
with the return = "archive"
argument) containing particle histories for connected particles or (b) a pf_path-class
object containing reconstructed paths must be supplied. Depending on the implementation, using a subset, all or an expanded sample of sampled locations, the functions apply KUD smoother(s) via a user-supplied estimation routine (i.e., kernelUD
or kud_around_coastline
). The functions extract the utilisation distribution(s) as raster
(s), combine distribution(s) (if necessary), apply a spatial mask (e.g. coastline), plot the processed distribution (if specified) and return this as a raster
.
pf_kud_1(
xpf,
bathy,
sample_size = NULL,
estimate_ud = adehabitatHR::kernelUD,
grid,
...,
scale = FALSE,
plot_by_time = FALSE,
prompt = TRUE,
chunks = 1L,
cl = NULL,
varlist = NULL,
mask = NULL,
plot = TRUE,
verbose = TRUE
)
pf_kud_2(
xpf,
bathy,
sample_size = NULL,
estimate_ud = adehabitatHR::kernelUD,
grid,
...,
mask = NULL,
plot = TRUE,
verbose = TRUE
)
xpf |
A |
bathy |
A |
sample_size |
(optional) An integer that defines the number of particles to sample from (a) particle histories or (b) each path in |
estimate_ud |
A function (either |
grid , ... |
Arguments passed to |
scale |
For |
plot_by_time , prompt |
For |
chunks , cl , varlist |
For |
mask |
(optional) A spatial mask (see |
plot |
A logical input that defines whether or not to plot the KUD. |
verbose |
A logical input that defines whether or not to print messages to the console to monitor function progress. |
These function create smooth KUD representations of particle or path samples from a PF algorithm (see pf_plot_map
). Two different methods are implemented.
pf_kud_1
implements KUD estimation by fitting a single kernel to sampled locations at each time step; time step-specific kernels are then summed and re-normalised. This method can be implemented for (a) the sampled particles that formed continuous paths from the start to the end of the time series (see pf_simplify
) or (b) reconstructed paths. The method has two main advantages. First, the kernel bandwidth is allowed to vary through time. Second, this implementation permits a Bayesian-style resampling process that can be used to account for particle uncertainty. Specifically, at each time step, large numbers of particles can be re-sampled, in line with their probability and with replacement, from the initial list of sampled particles; more probable locations are sampled more often and consequently have more influence on the KUD for each time step, thus accounting for particle probability.
A limitation with this method is that the fitting a KUD to each time step can be memory intensive and computationally intensive. To minimise memory requirements, by default (chunks = 1L
), the function starts with a blank map and iterates over each time step, sequentially adding KUDs to the map at each step. By continuously updating a single map, this option minimises memory requirements but is slow. A faster option is to split the time series into chunks, implement an iterative option within each chunk, and then join maps for each chunk. The advantage of this option is that memory use remains limited while computation time can be improved by the parallel processing of each chunk. This is implemented via chunks
, cl
and varlist
.
pf_kud_2
implements KUD estimation by fitting a single kernel to all sampled locations or to each path (and then aggregating KUDs across paths). The main advantage of this method is speed: unlike pf_kud_1
, KUDs are not fitted to the locations for each time step. The limitations are that kernel bandwidth is constant for all time steps and in most situations re-sampling cannot be used in the same way to account for particle uncertainty (due to memory limitations).
For particle-based implementations, this method is designed to be implemented for the subset of unique, sampled particles that formed continuous paths from the start to the end of the time series (see pf_simplify
and pf_plot_map
). These particles are used for KUD estimation. By default, all particles are used, but n =
sample_size
particles can be sampled at random, in line with their probability, if specified, for faster KUD estimation. Selected particles are then used to estimate a KUD by effectively treating samples as ‘relocations’, ultimately via kernelUD
. This distribution is then processed, plotted and returned.
For path-based implementations, this function is designed to be implemented for paths reconstructed by pf_simplify
. As for the particle-based implementation, for each path, all locations, or random sample of sample_size
locations are used to estimate a KUD by treating sampled locations as ‘relocations’. KUDs are processed and combined across paths into a single, average KUD. The advantage of this approach is that the overall probability of the paths can be incorporated in the estimation procedure via sampling or weights when path-specific KUDs are averaged (although that is not yet implemented).
The functions return a raster
of the KUD.
Edward Lavender
pf
, pf_simplify
, pf_plot_map
, kernelUD
, kud_around_coastline
, kud_around_coastline_fast
, eval_by_kud
#### Define particle samples for smoothing
# To do this, we will re-implement the pf() for the example dat_dcpf_histories
# ... dataset, but with a larger number of particles. This is necessary
# ... because kernel smoothing is only appropriate if there are
# ... enough locations to permit smoothing.
set.seed(1)
dcpf_args <- dat_dcpf_histories$args
dcpf_args$calc_distance_euclid_fast <- TRUE
dcpf_args$n <- 250L
out_dcpf_particles <- do.call(pf, dcpf_args)
#### Process particles and paths
out_dcpf_particles_1 <-
pf_simplify(out_dcpf_particles, return = "archive")
out_dcpf_particles_2 <-
pf_simplify(out_dcpf_particles, summarise_pr = TRUE, return = "archive")
out_dcpf_paths <-
pf_simplify(out_dcpf_particles, max_n_paths = 100L)
#### Define a grid across which to implement estimation
# This grid takes values of 0 on land and values of 1 in the sea
bathy <- out_dcpf_particles$args$bathy
grid <- raster::raster(raster::extent(bathy), nrows = 100, ncols = 100)
raster::values(grid) <- 0
grid <- raster::mask(grid, dat_coast, updatevalue = 1)
grid <- methods::as(grid, "SpatialPixelsDataFrame")
#### Example (1): Implement pf_kud_1() using default options
## Implementation based on particles
pp <- par(mfrow = c(1, 2), mar = c(3, 3, 3, 3))
pf_kud_1(out_dcpf_particles_1,
bathy = bathy,
sample_size = 500,
estimate_ud = kud_around_coastline_fast, grid = grid
)
## Implementation based on paths
if (flapper_run_parallel) {
pf_kud_1(out_dcpf_paths,
bathy = bathy,
sample_size = 500,
estimate_ud = kud_around_coastline_fast, grid = grid
)
prettyGraphics::add_sp_path(
x = out_dcpf_paths$cell_x,
y = out_dcpf_paths$cell_y,
length = 0.01
)
}
par(pp)
#### Example (2): Implement pf_kud_2() using default options
## Implementation based on particles
pp <- par(mfrow = c(1, 2), mar = c(3, 3, 3, 3))
pf_kud_2(out_dcpf_particles_2,
bathy = bathy,
estimate_ud = kud_around_coastline, grid = grid
)
## Implementation based on paths
pf_kud_2(out_dcpf_paths,
bathy = bathy,
estimate_ud = kud_around_coastline, grid = grid
)
prettyGraphics::add_sp_path(
x = out_dcpf_paths$cell_x, y = out_dcpf_paths$cell_y,
length = 0.01
)
par(pp)
#### Example (3): For improved speed with pf_kud_1(), use parallelisation
## Implementation based on particles
pp <- par(mfrow = c(1, 2), mar = c(3, 3, 3, 3))
pf_kud_1(out_dcpf_particles_1,
bathy = bathy,
sample_size = 500,
estimate_ud = flapper::kud_around_coastline_fast, grid = grid,
chunks = 2L,
cl = parallel::makeCluster(2L)
)
## Implementation based on paths
if (flapper_run_parallel) {
cl <- parallel::makeCluster(2L)
parallel::clusterEvalQ(cl = cl, library(raster))
pf_kud_1(out_dcpf_paths,
bathy = bathy,
sample_size = 500,
estimate_ud = flapper::kud_around_coastline_fast, grid = grid,
chunks = 2L,
cl = cl
)
}
par(pp)
#### Example (4): For improved speed with pf_kud_2(), use sample_size
## Implementation based on particles
pp <- par(mfrow = c(1, 2), mar = c(3, 3, 3, 3))
pf_kud_2(out_dcpf_particles_2,
bathy = bathy,
sample_size = 50, # sample 50 particles overall
estimate_ud = kud_around_coastline, grid = grid
)
## Implementation based on paths
pf_kud_2(out_dcpf_paths,
bathy = bathy,
sample_size = 50, # sample 50 particles per path
estimate_ud = kud_around_coastline, grid = grid
)
par(pp)
#### Example (5): Compare pf_kud_1() and pf_kud_2()
pp <- par(mfrow = c(2, 2), mar = c(3, 3, 3, 3))
kud_1a <- pf_kud_1(
xpf = out_dcpf_particles_1,
bathy = out_dcpf_particles$args$bathy,
sample_size = out_dcpf_particles$args$n,
estimate_ud = kud_around_coastline_fast, grid = grid,
plot = TRUE
)
kud_1b <- pf_kud_1(
xpf = out_dcpf_particles_1,
bathy = out_dcpf_particles$args$bathy,
sample_size = 500,
estimate_ud = kud_around_coastline_fast, grid = grid,
plot = TRUE
)
kud_1c <- pf_kud_1(
xpf = out_dcpf_particles_1,
bathy = out_dcpf_particles$args$bathy,
sample_size = 5000,
estimate_ud = kud_around_coastline_fast, grid = grid,
plot = TRUE
)
kud_1d <- pf_kud_2(
xpf = out_dcpf_particles_2,
bathy = out_dcpf_particles$args$bathy,
estimate_ud = kud_around_coastline, grid = grid,
plot = TRUE
)
par(pp)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.