| detector_create | R Documentation |
Creates an online (sequential) changepoint detector object that provides
a step-by-step interface to the FOCuS algorithm. Each call to
detector_update() adds new data, and get_statistics() computes
the current test statistic and detection result.
detector_create(
type,
dim_indexes = NULL,
quantiles = NULL,
pruning_mult = 2L,
pruning_offset = 1L,
side = "right",
anomaly_intensity = NULL,
rho = NULL,
mu0_arp = NULL
)
type |
Character string specifying detector type. One of:
|
dim_indexes |
List of integer vectors specifying projection index sets
for high-dimensional multivariate detectors. Not required for sequences of dimentions less than 5.
Each element is a vector of 0-based column indices. Default is |
quantiles |
Numeric vector of quantiles for nonparametric
( |
pruning_mult |
Integer. Candidate pruning multiplier parameter. Default is 2. |
pruning_offset |
Integer. Candidate pruning offset parameter. Default is 1. |
side |
Character string. For one-sided detectors, either |
anomaly_intensity |
Numeric scalar. Anomaly intensity threshold for
pruning candidates. Only candidates with sufficient signal magnitude are
retained. Default is |
rho |
Numeric vector. AR coefficients for AutoRegressive Process (ARP)
detectors. Required when |
mu0_arp |
Numeric scalar. Pre-change mean for ARP detectors (optional).
When provided, enables more efficient pruning by filtering candidates based on
the known pre-change parameter. Only used when |
The detector maintains sufficient statistics internally and uses pruning
to efficiently track candidate changepoints. The pruning_mult and
pruning_offset parameters control the pruning strategy.
AutoRegressive Process (ARP):
When type = "arp", the rho parameter must be provided as a numeric
vector of AR coefficients (lag-1, lag-2, ..., lag-p). The detector then computes
statistics optimal for detecting changepoints in AR(p) processes. Use
get_statistics(family = "arp") to retrieve the test statistics.
The optional theta0 parameter specifies the pre-change mean and is tied
to the pruning logic: if provided, it enables more efficient pruning by allowing
the algorithm to filter candidates based on the known pre-change parameter.
High-dimentional multivariate detectors:
For high-dimensional multivariate detection, computing the full hull would be too prohibitive.
Additionally, the complexity is expected to be log(n)^p, where n is the number of iterations and p is the
dimensions. So for short, high-dimentional sequences, it is possible to reconstruct the set of changepoint
locations by approximating the hull on projections in smaller dimentions. dim_indexes specifies which dimensions
to use for each projection of the convex hull for pruning. Use generate_projection_indexes() to
generate systematic projection sets.
NPFOCuS:
For non-parametric detection one needs to set detector(type = "npfocus") and the cost can be computed as get_statistics(family = "npfocus").
Any other cost will not work with this detector type. The quantiles vector argument is required, see get_statistics() for details.
An external pointer (SEXP) to the detector object. This should be
passed to other detector functions like detector_update() and
get_statistics().
# Univariate detector
det <- detector_create(type = "univariate")
detector_update(det, 0.5)
detector_update(det, 1.2)
r <- get_statistics(det, family = "gaussian")
print(r)
## Online (sequential) example
# Generate data with a changepoint
set.seed(123)
Y <- c(rnorm(500, mean = 0), rnorm(500, mean = 1))
det <- detector_create(type = "univariate")
stat_trace <- numeric(length(Y))
threshold <- 20
for (i in seq_along(Y)) {
detector_update(det, Y[i])
r <- get_statistics(det, family = "gaussian")
stat_trace[i] <- r$stat
if (!is.null(r$stat) && r$stat > threshold) {
cat("Online detection at", i, "estimate tau =", r$changepoint, "\n")
plot(stat_trace[1:i], type = "l", ylab = "Test Statistic", xlab = "Time")
break
}
}
# Multivariate detector with projections
dim_indexes <- list(c(0,1), c(1,2), c(0,2)) # 0-based indices
det_mv <- detector_create(type = "multivariate", dim_indexes = dim_indexes)
detector_update(det_mv, c(0.5, 1.2, -0.3))
# Nonparametric detector
quants <- qnorm(c(0.25, 0.5, 0.75))
det_np <- detector_create(type = "npfocus", quantiles = quants)
# One-sided univariate detector
det_one_sided <- detector_create(type = "univariate_one_sided", side = "left")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.