View source: R/check_model_fit.R
| check_model_fit | R Documentation |
Performs global envelope tests for nonparametric L, F, G, J, E, and V summary
functions (spatstat/GET).
These tests assess goodness-of-fit of the estimated model relative to a reference marked point pattern.
The reference marked point pattern can be supplied directly via reference_data (a marked ppp object),
or derived internally from a ldmppr_fit object.
check_model_fit(
reference_data = NULL,
t_min = 0,
t_max = 1,
process = c("self_correcting"),
process_fit = NULL,
anchor_point = NULL,
raster_list = NULL,
scaled_rasters = FALSE,
mark_model = NULL,
xy_bounds = NULL,
include_comp_inds = FALSE,
thinning = TRUE,
edge_correction = "none",
competition_radius = 15,
n_sim = 2500,
save_sims = TRUE,
verbose = TRUE,
seed = 0,
parallel = FALSE,
num_cores = max(1L, parallel::detectCores() - 1L),
set_future_plan = FALSE,
mark_mode = NULL,
fg_correction = c("km", "rs"),
max_attempts = NULL
)
reference_data |
(optional) a marked |
t_min |
minimum value for time. |
t_max |
maximum value for time. |
process |
type of process used (currently supports |
process_fit |
either an |
anchor_point |
(optional) vector of (x,y) coordinates of the point to condition on.
If |
raster_list |
(optional) list of raster objects used for mark prediction.
Required when |
scaled_rasters |
|
mark_model |
a mark model object used when |
xy_bounds |
(optional) vector of bounds as |
include_comp_inds |
|
thinning |
|
edge_correction |
type of edge correction to apply ( |
competition_radius |
positive numeric distance used when |
n_sim |
number of simulated datasets to generate. |
save_sims |
|
verbose |
|
seed |
integer seed for reproducibility. |
parallel |
|
num_cores |
number of workers to use when |
set_future_plan |
|
mark_mode |
(optional) mark generation mode: |
fg_correction |
correction used for F/G/J summaries ( |
max_attempts |
maximum number of simulation attempts when rejection occurs due to non-finite summaries. |
This function relies on the spatstat package for the calculation of the point pattern metrics
and the GET package for the global envelope tests. The L, F, G, J, E, and V functions are a collection of
non-parametric summary statistics that describe the spatial distribution of points and marks in a point pattern.
See the documentation for Lest(), Fest(), Gest(),
Jest(), Emark(), and Vmark() for more information.
Also, see the global_envelope_test() function for more information on the global envelope tests.
an object of class "ldmppr_model_check".
Baddeley, A., Rubak, E., & Turner, R. (2015). *Spatial Point Patterns: Methodology and Applications with R*. Chapman and Hall/CRC Press, London. ISBN 9781482210200. Available at: https://www.routledge.com/Spatial-Point-Patterns-Methodology-and-Applications-with-R/Baddeley-Rubak-Turner/p/book/9781482210200.
Myllymäki, M., & Mrkvička, T. (2023). GET: Global envelopes in R. arXiv:1911.06583 [stat.ME]. \Sexpr[results=rd]{tools:::Rd_expr_doi("10.48550/arXiv.1911.06583")}.
# Note: The example below is provided for illustrative purposes and may take some time to run.
data(small_example_data)
file_path <- system.file("extdata", "example_mark_model.rds", package = "ldmppr")
mark_model <- load_mark_model(file_path)
raster_paths <- list.files(system.file("extdata", package = "ldmppr"),
pattern = "\\.tif$", full.names = TRUE)
raster_paths <- raster_paths[!grepl("_med\\.tif$", raster_paths)]
rasters <- lapply(raster_paths, terra::rast)
scaled_raster_list <- scale_rasters(rasters)
reference_data <- generate_mpp(
locations = small_example_data[, c("x", "y")],
marks = small_example_data$size,
xy_bounds = c(0, 25, 0, 25)
)
estimated_parameters <- c(
0.05167978, 8.20702166, 0.02199940, 2.63236890,
1.82729512, 0.65330061, 0.86666748, 0.04681878
)
# Keep parallel=FALSE in examples to avoid setup overhead.
example_model_fit <- check_model_fit(
reference_data = reference_data,
t_min = 0,
t_max = 1,
process = "self_correcting",
process_fit = estimated_parameters,
raster_list = scaled_raster_list,
scaled_rasters = TRUE,
mark_model = mark_model,
xy_bounds = c(0, 25, 0, 25),
include_comp_inds = TRUE,
thinning = TRUE,
edge_correction = "none",
competition_radius = 10,
n_sim = 100,
save_sims = FALSE,
verbose = TRUE,
seed = 90210,
parallel = FALSE
)
plot(example_model_fit, which = 'combined')
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.