View source: R/sim-crw_in_polygon.r
crw_in_polygon | R Documentation |
Uses crw()
to simulate a random walk as series of equal-length steps
with turning angles drawn from a normal distribution inside a polygon.
crw_in_polygon(
polyg,
theta = c(0, 10),
stepLen = 100,
initPos = c(NA, NA),
initHeading = NA,
nsteps = 30,
inputCRS = NA,
cartesianCRS = NA,
sp_out = TRUE,
show_progress = TRUE
)
polyg |
A spatial polygon object of class |
theta |
A 2-element numeric vector with turn angle parameters ( |
stepLen |
A numeric scalar with total distance moved in each step, in meters. |
initPos |
A 2-element numeric vector with initial position
( |
initHeading |
A numeric scalar with initial heading in degrees. E.g., 0 = North; 90 = East, 180 = South, 270 = West; etc. |
nsteps |
A numeric scalar with number of steps to simulate. |
inputCRS |
A |
cartesianCRS |
Coordinate reference system used for simulations. Must be a Cartesian (projected) coordinate system. Must be given when input CRS is non-Cartesian (e.g., long-lat); optional otherwise. See Note. |
sp_out |
Logical. If TRUE (default) then output is an |
show_progress |
Logical. Progress bar and status messages will be shown if TRUE (default) and not shown if FALSE. |
If initPos = NA, then a starting point is randomly selected within
the polygon boundary. A path is simulated forward using crw()
.
Initial heading is also randomly selected if initHeading = NA
. When
a step crosses the polygon boundary, a new heading for that step is drawn
and the turn angle standard deviation is enlarged slightly for each
subsequent point that lands outside the polygon.
If input polyg
object is a data.frame with x and y columns and
sp_out = TRUE
, then output object
coordinate system is defined by inputCRS
. Coordinate system on output
will be same as input if polyg
contains a valid CRS.
When sp_out = TRUE
, an sf
object containing one
POINT
feature for each vertex in the simulated path.
OR
When sp_out = FALSE
, a two-column data frame containing:
x |
x coordinates |
y |
y coordinates |
in the same units as
polyg
.
The path is constructed in segments based on the minimum distance between the previous point and the closest polygon boundary.
Simulations are conducted within the coordinate system specified by
argument cartesianCRS
.
EPSG 3175 (cartesianCRS = 3175
) is recommended projected
coordinate system for the North American Great Lakes Basin and St. Lawrence
River system.
https://spatialreference.org/ref/epsg/nad83-great-lakes-and-st-lawrence-albers/.
C. Holbrook cholbrook@usgs.gov
crw, transmit_along_path, detect_transmissions
# Example 1 - data.frame input
mypolygon <- data.frame(x = c(-50, -50, 50, 50), y = c(-50, 50, 50, -50))
path_df <- crw_in_polygon(mypolygon,
theta = c(0, 20), stepLen = 10,
initPos = c(0, 0), initHeading = 0, nsteps = 50, sp_out = FALSE
)
class(path_df) # note object is data.frame
plot(path_df,
type = "o", pch = 20, asp = c(1, 1),
xlim = range(mypolygon$x), ylim = range(mypolygon$y)
)
polygon(mypolygon, border = "red")
# Example 2 - data.frame input; input CRS specified
mypolygon <- data.frame(
x = c(-84, -85, -85, -84),
y = c(45, 44, 45, 45)
)
path_df <- crw_in_polygon(mypolygon,
theta = c(0, 20),
stepLen = 1000,
initPos = c(-84.75, 44.75),
initHeading = 0,
nsteps = 50,
inputCRS = 4326,
cartesianCRS = 3175,
sp_out = FALSE
)
plot(path_df,
type = "o", pch = 20, asp = c(1, 1),
xlim = range(mypolygon$x), ylim = range(mypolygon$y)
)
class(path_df) # note object is data.frame
polygon(mypolygon, border = "red")
# Example 3 - sf POLYGON input
data(great_lakes_polygon)
# simulate in great lakes polygon
path_sf <- crw_in_polygon(great_lakes_polygon,
theta = c(0, 25),
stepLen = 10000,
initHeading = 0,
nsteps = 100,
cartesianCRS = 3175
)
# plot
plot(sf::st_geometry(great_lakes_polygon),
col = "lightgrey",
border = "grey"
)
points(sf::st_coordinates(path_sf), type = "o", pch = 20, col = "red")
# zoom in
plot(sf::st_geometry(great_lakes_polygon),
col = "lightgrey",
xlim = sf::st_bbox(path_sf)[c("xmin", "xmax")],
ylim = sf::st_bbox(path_sf)[c("ymin", "ymax")]
)
points(sf::st_coordinates(path_sf), type = "o", pch = 20, col = "red")
# Example 4 - SpatialPolygonsDataFrame input
data(greatLakesPoly)
# simulate in great lakes polygon
path_sp <- crw_in_polygon(greatLakesPoly,
theta = c(0, 25),
stepLen = 10000,
initHeading = 0,
nsteps = 100,
cartesianCRS = 3175,
sp_out = TRUE
)
# plot
plot(sf::st_as_sfc(greatLakesPoly), col = "lightgrey", border = "grey")
points(sf::st_coordinates(sf::st_as_sf(path_sp)),
type = "o", pch = 20,
col = "red"
)
# zoom in
plot(sf::st_as_sfc(greatLakesPoly),
col = "lightgrey", border = "grey",
xlim = sf::st_bbox(path_sp)[c("xmin", "xmax")],
ylim = sf::st_bbox(path_sp)[c("ymin", "ymax")]
)
points(sf::st_coordinates(sf::st_as_sf(path_sp)),
type = "o", pch = 20,
col = "red"
)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.