Title: Simulating Nonhomogeneous Poisson Point Processes
Version: 1.0.2
Description: Simulates events from one dimensional nonhomogeneous Poisson point processes (NHPPPs) as per Trikalinos and Sereda (2024, <doi:10.48550/arXiv.2402.00358> and 2024, <doi:10.1371/journal.pone.0311311>). Functions are based on three algorithms that provably sample from a target NHPPP: the time-transformation of a homogeneous Poisson process (of intensity one) via the inverse of the integrated intensity function (Cinlar E, "Theory of stochastic processes" (1975, ISBN:0486497996)); the generation of a Poisson number of order statistics from a fixed density function; and the thinning of a majorizing NHPPP via an acceptance-rejection scheme (Lewis PAW, Shedler, GS (1979) <doi:10.1002/nav.3800260304>).
License: GPL (≥ 3)
Imports: lifecycle, rstream, Rcpp (≥ 1.0.12)
LinkingTo: Rcpp
Encoding: UTF-8
Language: es
RoxygenNote: 7.3.2
Suggests: data.table, knitr, rlecuyer, rmarkdown, testthat, tictoc, truncnorm, withr
Config/Needs/website: rmarkdown
URL: https://bladder-ca.github.io/nhppp/, https://github.com/bladder-ca/nhppp
BugReports: https://github.com/bladder-ca/nhppp/issues
VignetteBuilder: knitr
Depends: R (≥ 2.10)
LazyData: true
NeedsCompilation: yes
Packaged: 2025-01-09 16:00:42 UTC; tom
Author: Thomas Trikalinos ORCID iD [aut, cre, cph], Yuliia Sereda ORCID iD [aut]
Maintainer: Thomas Trikalinos <thomas_trikalinos@brown.edu>
Repository: CRAN
Date/Publication: 2025-01-09 16:20:02 UTC

nhppp: Simulating Nonhomogeneous Poisson Point Processes

Description

logo

Simulates events from one dimensional nonhomogeneous Poisson point processes (NHPPPs) as per Trikalinos and Sereda (2024, doi:10.48550/arXiv.2402.00358 and 2024, doi:10.1371/journal.pone.0311311). Functions are based on three algorithms that provably sample from a target NHPPP: the time-transformation of a homogeneous Poisson process (of intensity one) via the inverse of the integrated intensity function (Cinlar E, "Theory of stochastic processes" (1975, ISBN:0486497996)); the generation of a Poisson number of order statistics from a fixed density function; and the thinning of a majorizing NHPPP via an acceptance-rejection scheme (Lewis PAW, Shedler, GS (1979) doi:10.1002/nav.3800260304).

Author(s)

Maintainer: Thomas Trikalinos thomas_trikalinos@brown.edu (ORCID) [copyright holder]

Authors:

See Also

Useful links:


Definite integral of l = exp(intercept + slope*t) at time t with L(t0) = 0

Description

Definite integral of l = exp(intercept + slope*t) starting at t0 – only for ⁠l+⁠.

Usage

Lambda_exp_form(t, intercept, slope, t0)

Arguments

t

(double) the time point

intercept

(double) the intercept

slope

(double) the slope

t0

(double) the starting time


Inverse of the definite integral of l = exp(intercept + slope*t) at time t

Description

Inverse of the definite integral of l = exp(intercept + slope*t) only for ⁠l+⁠.

Usage

Lambda_inv_exp_form(z, intercept, slope, t0)

Arguments

z

(double) the value of integrated rate for which you want to find the time

intercept

(double) the intercept

slope

(double) the slope

t0

(double) the starting time


Inverse of the definite integral of l = intercept + slope*t at time t

Description

Inverse of the definite integral of l = intercept + slope*t only for ⁠l+⁠.

Usage

Lambda_inv_linear_form(z, intercept, slope, t0)

Arguments

z

(double) the value of integrated rate for which you want to find the time

intercept

(double) the intercept

slope

(double) the slope

t0

(double) the starting time


Definite integral of l = intercept + slope*t at time t with L(t0) = 0

Description

Definite integral of l = intercept + slope*t starting at t0 – only for ⁠l+⁠.

Usage

Lambda_linear_form(t, intercept, slope, t0)

Arguments

t

(double) the time point

intercept

(double) the intercept

slope

(double) the slope

t0

(double) the starting time


Human mortality database age and sex specific rates for all cause deaths

Description

This is the 2015 annual death rates from the 2023 version of the Human Mortality Database for the USA.

Usage

annual_mortality_rates_2015

Format

annual_mortality_rates_2015

a data.table with 3 rows and 113 columns.

birth_cohort

Birth cohort as a calendar year

sex

sex: female, male, total.

age_0, ..., age 110+

age-specific death rates

Source

https://www.mortality.org


Check the validity of ppp samples arranged in matrix format

Description

Standard checks for a matrix of ordered times (event series in rows, times in columns). Check that the times in the columns are sorted, have unique values in ⁠[t_min, t_max]⁠, and has length size (if applicable).

Usage

check_ppp_sample_validity(
  times,
  t_min,
  t_max = NULL,
  size = NULL,
  atmost1 = FALSE,
  atleast1 = FALSE
)

Arguments

times

(vector, double | matrix) the times to be checked as vectors or matrices (time-vectors in rows)

t_min

(double | vector) the start of the time nterval

t_max

(double| vector) optional: the end of the time interval; if a vector, its length should match the number of rows of times.

size

(double) optional: the size of the vector

atmost1

(boolean) optional: at most one sample returned

atleast1

(boolean) optional: at least one sample returned

Value

None


Check the validity of a ppp vector.

Description

Standard checks for a vector of ordered times. Check that the times vector is sorted, has unique values, has all values in ⁠[t_min, t_max]⁠, and has length size (if applicable).

Usage

check_ppp_vector_validity(
  times,
  t_min,
  t_max = NULL,
  size = NULL,
  atmost1 = FALSE,
  atleast1 = FALSE
)

Arguments

times

(vector, double) the times to be checked

t_min

(double) the start of the time nterval

t_max

(double) optional: the end of the time interval

size

(double) optional: the size of the vector

atmost1

(boolean) optional: at most one sample returned

atleast1

(boolean) optional: at least one sample returned

Value

None


Check that two ppp vectors Q-Q agree

Description

Compare that the deciles of two vectors have absolute difference over average ratios less than threshold

Usage

compare_ppp_vectors(ppp1, ppp2, threshold = 0.15, showQQ = TRUE)

Arguments

ppp1

(vector, double) the first vector

ppp2

(vector, double) the second vector

threshold

(double) optional: the cutoff for a large absolute threshold

showQQ

(boolean) optional: show the QQ plot if the absolute value of the Difference vs Average ratio in any decile is bigger than the threshold

Value

None


Generic function for simulating from NHPPPs given the intensity function or the cumulative intensity function.

Description

This is a wrapper to the package's specific functions, and thus somewhat slower. For time-intensive simulations prefer one of the specific functions.

Usage

draw(
  Lambda = NULL,
  Lambda_inv = NULL,
  lambda = NULL,
  line_majorizer_intercept = NULL,
  line_majorizer_slope = NULL,
  line_majorizer_is_loglinear = FALSE,
  step_majorizer_vector = NULL,
  t_min = NULL,
  t_max = NULL,
  atmost1 = FALSE,
  atleast1 = FALSE
)

Arguments

Lambda

(function, double vector) the integrated (cumulative) rate of the NHPPP

Lambda_inv

(function, double vector) the inverse of ‘Lambda()’

lambda

(function) the instantaneous rate

line_majorizer_intercept

The intercept alpha of the loglinear majorizer function: alpha + beta * t or exp(alpha + beta * t)

line_majorizer_slope

The slope beta of the loglinear majorizer function: alpha + beta * t or exp(alpha + beta * t)

line_majorizer_is_loglinear

(boolean) if TRUE the majorizer is loglinear exp(alpha + beta * t); if FALSE it is a linear function

step_majorizer_vector

(vector, double) K constant majorizing rates, one per interval; all intervals are of equal length (regular)

t_min

(double) the lower bound of the interval

t_max

(double) the upper bound of the interval

atmost1

boolean, draw at most 1 event time

atleast1

boolean, draw at least 1 event time in interval

Value

a vector of event times


Simulate from a non homogeneous Poisson Point Process (NHPPP) over an interval when you know the cumulative intensity and its inverse.

Description

Sample NHPPP times using the inversion method

Usage

draw_cumulative_intensity(Lambda, Lambda_inv, t_min, t_max, atmost1 = FALSE)

Arguments

Lambda

(function, double vector) a continuous increasing R to R map which is the integrated rate of the NHPPP

Lambda_inv

(function, double vector) the inverse of Lambda()

t_min

(double) the lower bound of the time interval

t_max

(double) the upper bound of the time interval

atmost1

boolean, draw at most 1 event time

Value

a vector of event times (t_); if no events realize, a vector of length 0


Simulate from a non homogeneous Poisson Point Process (NHPPP) from (t_min, t_max) (inversion method)

Description

Sample NHPPP times using the inversion method, optionally using an rstream generator object

Usage

draw_cumulative_intensity_inversion(
  Lambda,
  Lambda_inv,
  t_min,
  t_max,
  atmost1 = FALSE
)

Arguments

Lambda

(function, double vector) a continuous increasing R to R map which is the integrated rate of the NHPPP

Lambda_inv

(function, double vector) the inverse of Lambda()

t_min

(double) the lower bound of the time interval

t_max

(double) the upper bound of the time interval

atmost1

boolean, draw at most 1 event time

Value

a vector of event times (t_); if no events realize, a vector of length 0


Simulate from a non homogeneous Poisson Point Process (NHPPP) from (t_min, t_max) (order statistics method)

Description

Sample NHPPP times using the order statistics method, optionally using an rstream generator object.

Usage

draw_cumulative_intensity_orderstats(
  Lambda,
  Lambda_inv,
  t_min,
  t_max,
  atmost1 = FALSE
)

Arguments

Lambda

(function, double vector) a continuous increasing R to R map which is the integrated rate of the NHPPP

Lambda_inv

(function, double vector) the inverse of Lambda()

t_min

(double) the lower bound of the time interval

t_max

(double) the upper bound of the time interval

atmost1

boolean, draw at most 1 event time

Value

a vector of event times (t_); if no events realize, a vector of length 0


Generic function for simulating from NHPPPs given the intensity function.

Description

Sample from NHPPPs given the intensity function This is a wrapper to the package's specific functions, and thus somewhat slower. For time-intensive simulations prefer one of the specific functions.

Usage

draw_intensity(
  lambda,
  line_majorizer_intercept = NULL,
  line_majorizer_slope = NULL,
  line_majorizer_is_loglinear = FALSE,
  step_majorizer_vector = NULL,
  t_min = NULL,
  t_max = NULL,
  atmost1 = FALSE
)

Arguments

lambda

(function) the instantaneous rate

line_majorizer_intercept

The intercept alpha of the loglinear majorizer function: alpha + beta * t or exp(alpha + beta * t)

line_majorizer_slope

The slope beta of the loglinear majorizer function: alpha + beta * t or exp(alpha + beta * t)

line_majorizer_is_loglinear

(boolean) if TRUE the majorizer is loglinear exp(alpha + beta * t); if FALSE it is a linear function

step_majorizer_vector

(vector, double) K constant majorizing rates, one per interval; all intervals are of equal length (regular)

t_min

(double) the lower bound of the interval

t_max

(double) the upper bound of the interval

atmost1

boolean, draw at most 1 event time

Value

a vector of event times


Simulate from a non homogeneous Poisson Point Process (NHPPP) from (t0, t_max) (thinning method)

Description

Sample NHPPP times using the thinning method

Usage

draw_intensity_line(
  lambda,
  majorizer_intercept,
  majorizer_slope,
  t_min,
  t_max,
  majorizer_is_loglinear = FALSE,
  atmost1 = FALSE
)

Arguments

lambda

(function) the instantaneous rate of the NHPPP.

majorizer_intercept

(double) the intercept (alpha) of the loglinear majorizer function.

majorizer_slope

(double) the slope (‘beta’) of the loglinear majorizer function.

t_min

(double) the lower bound of the time interval.

t_max

(double) the upper bound of the time interval.

majorizer_is_loglinear

(boolean) if TRUE the majorizer is loglinear exp(alpha + beta * t)

atmost1

boolean, draw at most 1 event time

Value

a vector of event times (t_); if no events realize, a vector of length 0


Simulate from a non homogeneous Poisson Point Process (NHPPP) from (t0, t_max) (thinning method) with piecewise constant_majorizer

Description

Sample NHPPP times using the thinning method

Usage

draw_intensity_step(lambda, majorizer_vector, time_breaks, atmost1 = FALSE)

Arguments

lambda

(function) the instantaneous rate of the NHPPP. A continuous function of time.

majorizer_vector

(scalar, double) K constant majorizing rates, one per interval

time_breaks

(vector, double) K+1 time points defining K intervals of constant rates: ⁠[t_1 = range_t[1], t_2)⁠: the first interval ⁠[t_k, t_{k+1})⁠: the k-th interval ⁠[t_{K}, t_{K+1} = range_t[2])⁠: the K-th (last) interval

atmost1

boolean, draw at most 1 event time

Value

a vector of event times (t_); if no events realize, a vector of length 0


Special case: Simulate from a non homogeneous Poisson Point Process (NHPPP) from (t_min, t_max) with linear intensity function (inversion method)

Description

Sample NHPPP times from a linear intensity function using the inversion method, optionally using an rstream generator

Usage

draw_sc_linear(intercept, slope, t_min, t_max, atmost1 = FALSE)

Arguments

intercept

(double) the intercept

slope

(double) the slope

t_min

(double) lower bound of the time interval

t_max

(double) upper bound of the time interval

atmost1

boolean, draw at most 1 event time

Value

a vector of event times (t_); if no events realize, a vector of length 0

Examples

x <- draw_sc_linear(intercept = 0, slope = 0.2, t_min = 0, t_max = 10)


Special case: Simulate from a non homogeneous Poisson Point Process (NHPPP) from (t_min, t_max) with log-linear intensity function (inversion method)

Description

Sample NHPPP times from an log linear intensity function using the inversion method, optionally using an rstream generator

Usage

draw_sc_loglinear(intercept, slope, t_min, t_max, atmost1 = FALSE)

Arguments

intercept

(double) the intercept in the exponent

slope

(double) the slope in the exponent

t_min

(double) lower bound of the time interval

t_max

(double) upper bound of the time interval

atmost1

boolean, draw at most 1 event time

Value

a vector of event times (t_); if no events realize, a vector of length 0

Examples

x <- draw_sc_loglinear(intercept = 0, slope = 0.2, t_min = 0, t_max = 10)


Simulate a piecewise constant-rate Poisson Point Process over ⁠(t_min, t_max]⁠ (inversion method) The intervals need not have the same length.

Description

Simulate a piecewise constant-rate Poisson Point Process over ⁠(t_min, t_max]⁠ (inversion method) The intervals need not have the same length.

Usage

draw_sc_step(lambda_vector, time_breaks, atmost1 = FALSE, atleast1 = FALSE)

Arguments

lambda_vector

(scalar, double) K constant rates, one per interval

time_breaks

(vector, double) K+1 time points defining K intervals of constant rates: ⁠[t_1 = range_t[1], t_2)⁠: the first interval ⁠[t_k, t_{k+1})⁠: the k-th interval ⁠[t_{K}, t_{K+1} = range_t[2])⁠: the K-th (last) interval

atmost1

boolean, draw at most 1 event time

atleast1

boolean, draw at least 1 event time

Value

a vector of event times t if no events realize, it will have 0 length

Examples

x <- draw_sc_step(lambda_vector = rep(1, 5), time_breaks = c(0:5))

Sampling from NHPPPs with piecewise constant intensities with same interval lengths (non-vectorized)

Description

Sampling from NHPPPs with piecewise constant intensities with same interval lengths (non-vectorized)

Usage

draw_sc_step_regular(
  Lambda_vector = NULL,
  lambda_vector = NULL,
  t_min = NULL,
  t_max = NULL,
  atmost1 = FALSE,
  atleast1 = FALSE
)

Arguments

Lambda_vector

(scalar, double) K integrated intensity rates at the end of each interval

lambda_vector

(scalar, double) K constant intensity rates, one per interval

t_min

(scalar, double) lower bound of the time interval

t_max

(scalar, double) upper bound of the time interval

atmost1

boolean, draw at most 1 event time

atleast1

boolean, draw at least 1 event time

Value

a vector of event times t if no events realize, it will have 0 length

Examples

x <- draw_sc_step_regular(Lambda_vector = 1:5, t_min = 0, t_max = 5)

Helper functions

Description

helper function that augments test_that::expect_no_error() to expect no error. Copied from the R6 source code.

Usage

expect_no_error(expr)

Arguments

expr

Expression.

Details

Small utility functions. Not to be exported to the user.


Piecewise constant (step) majorizer for K-Lipschitz functions over an interval (vectorized over the breaks argument).

Description

Return a piecewise constant (step) majorizer for K-Lipschitz functions over an interval. The function is vectorized over the breaks argument. The returned object has the same dimensions as breaks.

Usage

get_step_majorizer(fun, breaks, is_monotone = TRUE, K = 0)

Arguments

fun

A function object with a single argument x. If x is a matrix, fun should be vectorized over it.

breaks

(vector or matrix) The set of M+1 boundaries for the M subintervals in x. If breaks is a matrix, each row is treated as a separate set of breaks.

is_monotone

(boolean) Is the function monotone? (Default is TRUE.)

K

(double) A non-negative number for the Lipschitz cone. (Default is 0.)

Value

A vector of length M with the values of the piecewise constant majorizer

Examples

get_step_majorizer(fun = abs, breaks = -5:5, is_monotone = FALSE, K = 1)

Numerically evaluate the inverse of a function at a specific point

Description

Numerically evaluate the inverse of a function at a specific point

Usage

inverse_with_uniroot(
  f = f,
  y,
  min_x = 0,
  max_x = 1,
  min_y = f(min_x),
  max_y = f(max_x)
)

Arguments

f

(function) the function to be inverted; must be continuous and increasing

y

(scalar, double) the f(x)=y value in which to evaluate the inverse

min_x

(scalar, double) the min of the domain of f()

max_x

(scalar, double) the max of the domain of f()

min_y

(scalar, double) the min in the range of f()

max_y

(scalar, double) the max in the range of f()

Value

(scalar, double) vector of x=f^(-1)(y): the inverted value

Examples

inverse_with_uniroot(f = function(x) {
  2 * x
}, y = 0.5)

Numerically evaluate the inverse of a monotonically increasing continuous function from R to R at specific points.

Description

Numerically evaluate the inverse of a monotonically increasing continuous function from R to R at specific points.

Usage

inverse_with_uniroot_sorted(
  f,
  y,
  range_x = c(0, 10),
  range_y = c(f(range_x[1]), f(range_x[2]))
)

Arguments

f

(function) the function to be inverted; must be continuous and increasing

y

(vector, double) the f(x)=y values in which to evaluate the inverse; must be in ascending order

range_x

(vector, double) the min and max of the domain of f()

range_y

(vector, double) the min and max in the range of f()

Value

(vector, double) vector of x=f^(-1)(y): the inverted values

Examples

inverse_with_uniroot_sorted(f = function(x) {
  2 * x
}, y = c(0, 0.5))

Helper function for the vectorized versions of sampling functions. Takes the usual ways that lambda_mat and Lambda_mat are specified and returns Lambda_mat.

Description

Helper function for the vectorized versions of sampling functions. Takes the usual ways that lambda_mat and Lambda_mat are specified and returns Lambda_mat.

Usage

make_cumulative_Lambda_matrix(
  Lambda_mat = NULL,
  lambda_mat = NULL,
  interval_duration = NULL
)

Arguments

Lambda_mat

a matrix of cumulative intensities or missing

lambda_mat

a matrix of intensities or missing

interval_duration

a vector with the same number of elements as the rows of Lambda_mat

Value

A matrix (r x 2), row-expanded if needed


Helper function for the vectorized versions of sampling functions. Takes the usual ways that lambda_mat and Lambda_mat are specified and returns lambda_mat.

Description

Helper function for the vectorized versions of sampling functions. Takes the usual ways that lambda_mat and Lambda_mat are specified and returns lambda_mat.

Usage

make_lambda_matrix(
  lambda_mat = NULL,
  Lambda_mat = NULL,
  interval_duration = NULL
)

Arguments

lambda_mat

a matrix of intensities or missing

Lambda_mat

a matrix of cumulative intensities or missing

interval_duration

a vector with the same number of elements as the rows of Lambda_mat

Value

A matrix (r x 2), row-expanded if needed


Helper function for the vectorized versions of sampling functions. Takes the usual ways that range_t is specified (a 2-vector, a 1 x 2 or an r x 2 matrix) and returns a r x 2 matrix.

Description

Helper function for the vectorized versions of sampling functions. Takes the usual ways that range_t is specified (a 2-vector, a 1 x 2 or an r x 2 matrix) and returns a r x 2 matrix.

Usage

make_range_t_matrix(range_t, n_rows)

Arguments

range_t

a 2-vector, a 1 x 2 or an r x 2 matrix

n_rows

the number of rows in the fully expanded matrix (r)

Value

A matrix (r x 2), row-expanded if needed


Return matrix with column-wise cumulative sum No checks for arguments is done.

Description

Return matrix with column-wise cumulative sum No checks for arguments is done.

Usage

mat_cumsum_columns(X)

Arguments

X

(matrix)

Value

matrix


Return matrix with column-wise cumulative sum replacing cells larger than ceil with NA. No checks for arguments is done.

Description

Return matrix with column-wise cumulative sum replacing cells larger than ceil with NA. No checks for arguments is done.

Usage

mat_cumsum_columns_with_scalar_ceiling(X, ceil = Inf)

Arguments

X

(matrix)

ceil

(double or Inf) the ceiling to be applied

Value

matrix


Return matrix with column-wise cumulative sum replacing cells larger than ceil with NA. No checks for arguments is done.

Description

Return matrix with column-wise cumulative sum replacing cells larger than ceil with NA. No checks for arguments is done.

Usage

mat_cumsum_columns_with_vector_ceiling(X, ceil = Inf)

Arguments

X

(matrix)

ceil

(vector or Inf) the set of ceilings to be applied, per row of X

Value

matrix


Return matrix with column-wise differencing. No checks for arguments is done.

Description

Return matrix with column-wise differencing. No checks for arguments is done.

Usage

mat_diff_columns(X)

Arguments

X

(matrix)

Value

matrix


Usage: matrix_cumsum_columns( X )

Description

Usage: matrix_cumsum_columns( X )

Usage

matrix_cumsum_columns(X)

Usage: matrix_cumsum_columns_inplace( X )

Description

Usage: matrix_cumsum_columns_inplace( X )

Usage

matrix_cumsum_columns_inplace(X)

Usage: matrix_diff_columns( X )

Description

Usage: matrix_diff_columns( X )

Usage

matrix_diff_columns(X)

Usage: matrix_diff_columns_inplace( X )

Description

Usage: matrix_diff_columns_inplace( X )

Usage

matrix_diff_columns_inplace(X)

Simulate a homogeneous Poisson Point Process in (t_min, t_max]

Description

Simulate a homogeneous Poisson Point Process in (t_min, t_max]

Usage

ppp(rate, t_min, t_max, atmost1 = FALSE, tol = 10^-6)

Arguments

rate

(scalar, double) constant instantaneous rate

t_min

(scalar, double) the lower bound of the time interval

t_max

(scalar, double) the upper bound of the time interval

atmost1

boolean, draw at most 1 event time

tol

the probability that we will have more than the drawn events in (t_min, t_max]

Value

a vector of event times t if no events realize, it will have 0 length

Examples

x <- ppp(rate = 1, t_min = 0, t_max = 10, tol = 10^-6)

Simulate a homogeneous Poisson Point Process over (t_min, t_max] (order statistics method)

Description

Internal function – not to be exported. Same as ppp but uses the Order Statistics algorithm.

Usage

ppp2(rate, t_min, t_max, atmost1 = FALSE)

Arguments

rate

(scalar, double) constant instantaneous rate

t_min

(scalar, double) the lower bound of the time interval

t_max

(scalar, double) the upper bound of the time interval

atmost1

boolean, draw at most 1 event time

Value

a vector of event times t if no events realize, it will have 0 length

Examples

x <- ppp(rate = 1, t_min = 0, t_max = 10, tol = 10^-6)

Simulate exactly n points from a homogeneous Poisson Point Process over (t_min, t_max]

Description

Simulate exactly n points from a homogeneous Poisson Point Process over (t_min, t_max]

Usage

ppp_exactly_n(n, t_min, t_max)

Arguments

n

(int) the number of points to be simulated

t_min

(double) the lower bound of the time interval

t_max

(double) the upper bound of the time interval

Value

a vector of event times of size n

Examples

x <- ppp_exactly_n(n = 10, t_min = 0, t_max = 10)

Simulate specific number of points from a homogeneous Poisson Point Process over (t_min, t_max]

Description

[Deprecated] Use ppp_exactly_n instead.

Usage

ppp_n(size, range_t = c(0, 10), rng_stream = NULL)

Arguments

size

(int) the number of points to be simulated

range_t

(vector, double) min and max of the time interval

rng_stream

an rstream object

Value

a vector of event times of size size

Examples

x <- ppp_n(size = 10, range_t = c(0, 10))

Simulate n events from a homogeneous Poisson Point Process.

Description

Simulate n events from a homogeneous Poisson Point Process.

Usage

ppp_next_n(n = 1, rate = 1, t_min = 0, rng_stream = deprecated())

Arguments

n

scalar number of samples

rate

scalar instantaneous rate

t_min

scalar for the starting time value

rng_stream

[Deprecated] an rstream object

Value

a vector with event times t (starting from t_min)

Examples

x <- ppp_next_n(n = 10, rate = 1, t_min = 0)

Simulate a homogeneous Poisson Point Process over (t_min, t_max] (order statistics method)

Description

[Deprecated] Use ppp2 instead.

Usage

ppp_orderstat(range_t = c(0, 10), rate = 1, rng_stream = NULL, atmost1 = FALSE)

Arguments

range_t

(vector, double) min and max of the time interval

rate

(scalar, double) constant instantaneous rate

rng_stream

an rstream object

atmost1

boolean, draw at most 1 event time

Value

a vector of event times t if no events realize, it will have 0 length

Examples

x <- ppp_orderstat(range_t = c(0, 10), rate = 1)

Simulate a homogeneous Poisson Point Process over (t_min, t_max]

Description

[Deprecated] Use ppp instead.

Usage

ppp_sequential(
  range_t = c(0, 10),
  rate = 1,
  tol = 10^-6,
  rng_stream = NULL,
  atmost1 = FALSE
)

Arguments

range_t

(vector, double) min and max of the time interval

rate

(scalar, double) constant instantaneous rate

tol

the probability that we will have more than the drawn events in (t_min, t_max]

rng_stream

an rstream object

atmost1

boolean, draw at most 1 event time

Value

a vector of event times t if no events realize, it will have 0 length

Examples

x <- ppp_sequential(range_t = c(0, 10), rate = 1, tol = 10^-6)

Read code from text file as string

Description

Read code from text file as string

Usage

read_code(codeFile)

Arguments

codeFile

Path to file

Value

codeFile contents as a character string


Exponential random samples from rstream objects

Description

Sample from rstream objects

Usage

rng_stream_rexp(size = 1, rate = 1, rng_stream = NULL)

Arguments

size

Integer, number of samples

rate

Positive number, the rate (i.e., 1/mean)

rng_stream

(rstream) an rstream object or NULL

Value

a vector of exponential variates of size size

Examples

rng_stream_rexp(10)

Poisson random samples from rstream objects

Description

Sample from rstream objects

Usage

rng_stream_rpois(size = 1, lambda = 1, rng_stream = NULL)

Arguments

size

Integer, number of samples

lambda

Positive number, the mean

rng_stream

(rstream) an rstream object or NULL

Value

a vector of counts of size size

Examples

rng_stream_rpois(10)

Uniform random samples from rstream objects

Description

Sample from rstream objects

Usage

rng_stream_runif(size = 1, minimum = 0, maximum = 1, rng_stream = NULL)

Arguments

size

Integer, number of samples

minimum

Lower bound

maximum

Upper bound

rng_stream

(rstream) an rstream object or NULL

Value

a vector of uniform variates of size size

Examples

rng_stream_runif(10)

Zero-truncated Poisson random samples from rstream objects

Description

Sample from rstream objects

Usage

rng_stream_rztpois(size = 1, lambda = 1, rng_stream = NULL)

Arguments

size

Integer, number of samples

lambda

Positive number, the mean of the original (untruncated) Poisson distribution

rng_stream

(rstream) an rstream object or NULL

Value

a vector of non zero counts of size size

Examples

rng_stream_rztpois(10)

Zero-truncated Poisson random samples (basic R)

Description

Sample zero-truncated Poisson random samples (basic R)

Usage

rztpois(n, lambda)

Arguments

n

Integer, number of samples

lambda

Positive number, the mean of the original (untruncated) Poisson distribution

Value

a vector of non zero counts of size n

Examples

rztpois(10, 1)
rztpois(10, 1:10)

Simpson's method to integrate a univariate function.

Description

Simpson's method to integrate a univariate continuous function. Faster that R's integrate() and precise enough, but does not do any checks. The error is at most ⁠M (b-a)^5/(180 n^4)⁠ where M is the maximum of the fourth derivative of the integrand in the interval ⁠[a, b]⁠.

Usage

simpson_num_integr(f, a, b, n)

Arguments

f

function that takes a single argument

a

the lower limit of integration

b

the upper limit of integration

n

integer, number of integration points with a and b

Value

numeric, the integration value examples #expect 1 simpson_num_integr(sin, 0, pi/2, 100) #max error for simpson_num_integr(sin, 0, pi/2, 100) is 5.312842e-10 1 * (pi/2 - 0)^5/(180 * 100^4)


Vectorized generic function for simulating from NHPPPs given the intensity function or the cumulative intensity function

Description

This is a wrapper to the package's specific functions, and thus slightly slower. For time-intensive simulations prefer one of the specific functions.

Usage

vdraw(
  lambda = NULL,
  lambda_args = NULL,
  Lambda_maj_matrix = NULL,
  lambda_maj_matrix = NULL,
  Lambda = NULL,
  Lambda_inv = NULL,
  Lambda_args = NULL,
  Lambda_inv_args = NULL,
  t_min = NULL,
  t_max = NULL,
  rate_matrix_t_min = NULL,
  rate_matrix_t_max = NULL,
  tol = 10^-6,
  atmost1 = FALSE,
  atleast1 = FALSE,
  atmostB = NULL
)

Arguments

lambda

(function) intensity function, vectorized

lambda_args

(list) optional arguments to pass to lambda

Lambda_maj_matrix

(matrix) integrated intensity rates at the end of each interval

lambda_maj_matrix

(matrix) intensity rates, one per interval

Lambda

(function, double vector) an increasing function which is the integrated rate of the NHPPP. It should take a vectorized argument t for times and an optional arguments list.

Lambda_inv

(function, double vector) the inverse of Lambda(), also in vectorized form It should take a vectorized argument z and an optional arguments list.

Lambda_args

(list) optional arguments to pass to Lambda.

Lambda_inv_args

(list) optional arguments to pass to Lambda_inv().

t_min

(scalar | vector | column matrix) is the lower bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_min.

t_max

(scalar | vector | column matrix) is the upper bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_max.

rate_matrix_t_min

(scalar | vector | column matrix) is the lower bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

rate_matrix_t_max

(scalar | vector | column matrix) the upper bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

tol

(scalar, double) tolerance for the number of events

atmost1

boolean, draw at most 1 event time

atleast1

boolean, draw at least 1 event time

atmostB

If not NULL, draw at most B (B>0) event times. NULL means ignore.

Value

a vector of event times


Vectorized simulation from a non homogeneous Poisson Point Process (NHPPP) from (t_min, t_max) given the cumulative intensity function and its inverse

Description

Sample NHPPP times using the cumulative intensity function and its inverse.

Usage

vdraw_cumulative_intensity(
  Lambda,
  Lambda_inv,
  t_min,
  t_max,
  Lambda_args = NULL,
  Lambda_inv_args = NULL,
  tol = 10^-6,
  atmost1 = FALSE,
  atleast1 = FALSE
)

Arguments

Lambda

(function, double vector) an increasing function which is the integrated rate of the NHPPP. It should take a vectorized argument t for times and an optional arguments list.

Lambda_inv

(function, double vector) the inverse of Lambda(), also in vectorized form It should take a vectorized argument z and an optional arguments list.

t_min

(scalar | vector | column matrix) the lower bound of the interval for each sampled point process The length of this argument is the number of point processes that should be drawn.

t_max

(scalar | vector | column matrix) the upper bound of the interval for each sampled point process The length of this argument is the number of point processes that should be drawn.

Lambda_args

(list) optional arguments to pass to Lambda.

Lambda_inv_args

(list) optional arguments to pass to Lambda_inv().

tol

the tolerange for the calulations.

atmost1

boolean, draw at most 1 event time per sampled point process.

atleast1

boolean, draw at least 1 event time

Value

a matrix of event times with one row per sampled point process.


Vectorized sampling from a non homogeneous Poisson Point Process (NHPPP) from an interval (thinning method) with piecewise constant_majorizers (C++)

Description

Vectorized sampling from a non homogeneous Poisson Point Process (NHPPP) from an interval (thinning method) with piecewise constant_majorizers. The majorizers are step functions over equal-length time intevals.

Usage

vdraw_intensity(
  lambda = NULL,
  lambda_args = NULL,
  Lambda_maj_matrix = NULL,
  lambda_maj_matrix = NULL,
  rate_matrix_t_min = NULL,
  rate_matrix_t_max = NULL,
  t_min = NULL,
  t_max = NULL,
  tol = 10^-6,
  atmost1 = FALSE,
  atleast1 = FALSE,
  atmostB = NULL
)

Arguments

lambda

(function) intensity function, vectorized

lambda_args

(list) optional arguments to pass to lambda

Lambda_maj_matrix

(matrix) integrated intensity rates at the end of each interval

lambda_maj_matrix

(matrix) intensity rates, one per interval

rate_matrix_t_min

(scalar | vector | column matrix) is the lower bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

rate_matrix_t_max

(scalar | vector | column matrix) the upper bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

t_min

(scalar | vector | column matrix) is the lower bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_min.

t_max

(scalar | vector | column matrix) is the upper bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_max.

tol

(scalar, double) tolerance for the number of events

atmost1

boolean, draw at most 1 event time

atleast1

boolean, draw at least 1 event time

atmostB

If not NULL, draw at most B (B>0) event times. NULL means ignore.

Value

a matrix of event times (columns) per draw (rows) NAs are structural empty spots

Examples

x <- vdraw_intensity(
  lambda = function(x, ...) 0.1 * x,
  lambda_maj_matrix = matrix(rep(1, 5), nrow = 1),
  rate_matrix_t_min = 1,
  rate_matrix_t_max = 5
)

Vectorized sampling from a non homogeneous Poisson Point Process (NHPPP) given the intensity function with piecewise constant_majorizers (C++)

Description

Vectorized sampling from a non homogeneous Poisson Point Process (NHPPP) using the intensity function and piecewise constant_majorizers. The majorizers are step functions over equal-length time intevals.

Usage

vdraw_intensity_step_regular_cpp(
  lambda = NULL,
  lambda_args = NULL,
  Lambda_maj_matrix = NULL,
  lambda_maj_matrix = NULL,
  rate_matrix_t_min = NULL,
  rate_matrix_t_max = NULL,
  t_min = NULL,
  t_max = NULL,
  tol = 10^-6,
  atmost1 = FALSE,
  atmostB = NULL
)

Arguments

lambda

(function) intensity function, vectorized

lambda_args

(list) optional arguments to pass to lambda

Lambda_maj_matrix

(matrix) integrated intensity rates at the end of each interval

lambda_maj_matrix

(matrix) intensity rates, one per interval

rate_matrix_t_min

(scalar | vector | column matrix) is the lower bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

rate_matrix_t_max

(scalar | vector | column matrix) the upper bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

t_min

(scalar | vector | column matrix) is the lower bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_min.

t_max

(scalar | vector | column matrix) is the upper bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_max.

tol

(scalar, double) tolerance for the number of events

atmost1

boolean, draw at most 1 event time

atmostB

If not NULL, draw at most B (B>0) event times. NULL means ignore.

Value

a matrix of event times (columns) per draw (rows) NAs are structural empty spots


Vectorized sampling from a non homogeneous Poisson Point Process (NHPPP) from an interval (thinning method) with piecewise constant_majorizers (R) – but can be forced to sample from zero-truncated proposals.

Description

Vectorized sampling from a non homogeneous Poisson Point Process (NHPPP) from an interval (thinning method) with piecewise constant_majorizers. The majorizers are step functions over equal-length time intevals. This function is used for obtainning proposals for vztdraw_intensity_step_regular()

Usage

vdraw_intensity_step_regular_forcezt(
  lambda = NULL,
  lambda_args = NULL,
  Lambda_maj_matrix = NULL,
  lambda_maj_matrix = NULL,
  rate_matrix_t_min = NULL,
  rate_matrix_t_max = NULL,
  t_min = NULL,
  t_max = NULL,
  tol = 10^-6,
  atmost1 = FALSE,
  force_zt_majorizer = FALSE,
  ...
)

Arguments

lambda

(function) intensity function, vectorized

lambda_args

(list) optional arguments to pass to lambda

Lambda_maj_matrix

(matrix) integrated intensity rates at the end of each interval

lambda_maj_matrix

(matrix) intensity rates, one per interval

rate_matrix_t_min

(scalar | vector | column matrix) is the lower bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

rate_matrix_t_max

(scalar | vector | column matrix) the upper bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

t_min

(scalar | vector | column matrix) is the lower bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_min.

t_max

(scalar | vector | column matrix) is the upper bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_max.

tol

(scalar, double) tolerance for the number of events

atmost1

boolean, draw at most 1 event time

force_zt_majorizer

boolean, force the majorizer to be zero-truncated. This option is used when the function is called to make proposals for vztdraw_intensity_step_regular(). In general, do not set this option to TRUE.


Vectorized sampling from NHPPPs with piecewise constant intensities with same interval lengths

Description

Simulate a piecewise constant-rate Poisson Point Process over ⁠(t_min, t_max]⁠ (inversion method) where the intervals have the same length (are "regular").

Usage

vdraw_sc_step_regular(
  lambda_matrix = NULL,
  Lambda_matrix = NULL,
  rate_matrix_t_min = NULL,
  rate_matrix_t_max = NULL,
  t_min = NULL,
  t_max = NULL,
  tol = 10^-6,
  atmost1 = FALSE,
  atmostB = NULL,
  atleast1 = FALSE
)

Arguments

lambda_matrix

(matrix) intensity rates, one per interval

Lambda_matrix

(matrix) integrated intensity rates at the end of each interval

rate_matrix_t_min

(scalar | vector | column matrix) is the lower bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

rate_matrix_t_max

(scalar | vector | column matrix) the upper bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

t_min

(scalar | vector | column matrix) is the lower bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_min.

t_max

(scalar | vector | column matrix) is the upper bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_max.

tol

(scalar, double) tolerance for the number of events

atmost1

boolean, draw at most 1 event time

atmostB

If not NULL, draw at most B (B>0) event times. NULL means ignore.

atleast1

boolean, draw at least 1 event time

Value

a vector of event times t if no events realize, it will have 0 length

Examples

x <- vdraw_sc_step_regular(
  Lambda_matrix = matrix(1:5, nrow = 1),
  rate_matrix_t_min = 100,
  rate_matrix_t_max = 110,
  atmost1 = FALSE
)

Vectorized sampling from NHPPPs with piecewise constant intensities with same interval lengths (C++)

Description

Simulate a piecewise constant-rate Poisson Point Process over ⁠(t_min, t_max]⁠ (inversion method) where the intervals have the same length (are "regular").

Usage

vdraw_sc_step_regular_cpp(
  lambda_matrix = NULL,
  Lambda_matrix = NULL,
  rate_matrix_t_min = NULL,
  rate_matrix_t_max = NULL,
  t_min = NULL,
  t_max = NULL,
  tol = 10^-6,
  atmost1 = FALSE,
  atmostB = NULL
)

Arguments

lambda_matrix

(matrix) intensity rates, one per interval

Lambda_matrix

(matrix) integrated intensity rates at the end of each interval

rate_matrix_t_min

(scalar | vector | column matrix) is the lower bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

rate_matrix_t_max

(scalar | vector | column matrix) the upper bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

t_min

(scalar | vector | column matrix) is the lower bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_min.

t_max

(scalar | vector | column matrix) is the upper bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_max.

tol

(scalar, double) tolerance for the number of events

atmost1

boolean, draw at most 1 event time

atmostB

If not NULL, draw at most B (B>0) event times. NULL means ignore.

Value

a matrix of event times t, with rows corresponding to the sampled point processes.

Examples

x <- vdraw_sc_step_regular_cpp(
  Lambda_matrix = matrix(1:5, nrow = 1),
  rate_matrix_t_min = 100,
  rate_matrix_t_max = 110,
  atmost1 = FALSE
)

Vectorized sampling from a zero-truncated non homogeneous Poisson Point Process (NHPPP) from an interval (thinning method) with piecewise constant_majorizers

Description

Vectorized sampling from a zero-truncated non homogeneous Poisson Point Process (NHPPP) from an interval (thinning method) with piecewise constant_majorizers. The majorizers are step functions over equal-length time intevals.

Usage

vztdraw_intensity(
  lambda = NULL,
  lambda_args = NULL,
  Lambda_maj_matrix = NULL,
  lambda_maj_matrix = NULL,
  range_t = NULL,
  tol = 10^-6,
  atmost1 = FALSE,
  ...
)

Arguments

lambda

(function) a vectorized intensity function, with one or two arguments. The first is time. The optional second is a named list with additional arguments.

lambda_args

(list) optional list of named arguments for lambda()

Lambda_maj_matrix

(matrix) for the majorizeintegrated intensity rates at the end of each interval

lambda_maj_matrix

(matrix) intensity rates, one per interval

range_t

(vector, or matrix) t_min and t_max, possibly vectorized

tol

(scalar, double) tolerance for the number of events

atmost1

boolean, draw at most 1 event time

...

(any) other arguments (ignored – used for flexibility in calling from other functions)

Value

a matrix of event times (columns) per draw (rows) NAs are structural empty spots


Vectorized sampling from a zero-truncated non homogeneous Poisson Point Process (NHPPP) from an interval (thinning method) with piecewise constant_majorizers (R)

Description

Vectorized sampling from a zero-truncated non homogeneous Poisson Point Process (NHPPP) from an interval (thinning method) with piecewise constant_majorizers. The majorizers are step functions over equal-length time intevals.

Usage

vztdraw_intensity_step_regular(
  lambda = NULL,
  lambda_args = NULL,
  Lambda_maj_matrix = NULL,
  lambda_maj_matrix = NULL,
  rate_matrix_t_min = NULL,
  rate_matrix_t_max = NULL,
  t_min = NULL,
  t_max = NULL,
  atmost1 = FALSE,
  ...
)

Arguments

lambda

(function) intensity function, vectorized

lambda_args

(list) optional arguments to pass to lambda

Lambda_maj_matrix

(matrix) integrated intensity rates at the end of each interval

lambda_maj_matrix

(matrix) intensity rates, one per interval

rate_matrix_t_min

(scalar | vector | column matrix) is the lower bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

rate_matrix_t_max

(scalar | vector | column matrix) the upper bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

t_min

(scalar | vector | column matrix) is the lower bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_min.

t_max

(scalar | vector | column matrix) is the upper bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_max.

atmost1

boolean, draw at most 1 event time


Vectorized sampling from zero-truncated NHPPPs with piecewise constant intensities with same interval lengths

Description

Simulate a piecewise constant-rate Poisson Point Process over ⁠(t_min, t_max]⁠ (inversion method) where the intervals have the same length (are "regular").

Usage

vztdraw_sc_step_regular_cpp(
  lambda_matrix = NULL,
  Lambda_matrix = NULL,
  rate_matrix_t_min = NULL,
  rate_matrix_t_max = NULL,
  t_min = NULL,
  t_max = NULL,
  atmost1 = FALSE
)

Arguments

lambda_matrix

(matrix) intensity rates, one per interval

Lambda_matrix

(matrix) integrated intensity rates at the end of each interval

rate_matrix_t_min

(scalar | vector | column matrix) is the lower bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

rate_matrix_t_max

(scalar | vector | column matrix) the upper bound of the time interval for each row of (Lambda|lambda)_maj_matrix. The length of this argument is the number of point processes that should be drawn.

t_min

(scalar | vector | column matrix) is the lower bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_min.

t_max

(scalar | vector | column matrix) is the upper bound of a subinterval of (rate_matrix_t_min, rate_matrix_t_max]. If set, times are sampled from the subinterval. If omitted, it is equivalent to rate_matrix_t_max.

atmost1

boolean, draw at most 1 event time

Value

a vector of event times t if no events realize, it will have 0 length


Simulate from a zero-truncated non homogeneous Poisson Point Process (zt-NHPPP) from (t_min, t_max) (order statistics method)

Description

Sample zero-truncated NHPPP times using the order statistics method, optionally using an rstream generator

Usage

ztdraw_cumulative_intensity(Lambda, Lambda_inv, t_min, t_max, atmost1 = FALSE)

Arguments

Lambda

(function, double vector) a continuous increasing R to R map which is the integrated rate of the NHPPP

Lambda_inv

(function, double vector) the inverse of Lambda()

t_min

(double) the lower bound of the time interval

t_max

(double) the upper bound of the time interval

atmost1

(boolean) draw at most 1 event time

Value

a vector of at least 1 event times


Generic function for simulating from zero-truncated NHPPPs given the intensity function.

Description

Sample from zero-truncated NHPPP given the intensity function This is a wrapper to the package's specific functions, and thus somewhat slower. For time-intensive simulations prefer one of the specific functions.

Usage

ztdraw_intensity(
  lambda,
  line_majorizer_intercept = NULL,
  line_majorizer_slope = NULL,
  line_majorizer_is_loglinear = FALSE,
  step_majorizer_vector = NULL,
  t_min = NULL,
  t_max = NULL,
  atmost1 = FALSE
)

Arguments

lambda

(function) the instantaneous rate

line_majorizer_intercept

The intercept alpha of the loglinear majorizer function: alpha + beta * t or exp(alpha + beta * t)

line_majorizer_slope

The slope beta of the loglinear majorizer function: alpha + beta * t or exp(alpha + beta * t)

line_majorizer_is_loglinear

(boolean) if TRUE the majorizer is loglinear exp(alpha + beta * t); if FALSE it is a linear function

step_majorizer_vector

(vector, double) K constant majorizing rates, one per interval; all intervals are of equal length (regular)

t_min

(double) the lower bound of the interval

t_max

(double) the upper bound of the interval

atmost1

boolean, draw at most 1 event time

Value

a vector of at least 1 event times


Simulate size samples from a zero-truncated non homogeneous Poisson Point Process (zt-NHPPP) from (t0, t_max) (thinning method)

Description

Sample zero-truncated NHPPP intensity times using the thinning method

Usage

ztdraw_intensity_line(
  lambda,
  majorizer_intercept,
  majorizer_slope,
  t_min,
  t_max,
  majorizer_is_loglinear = FALSE,
  atmost1 = FALSE
)

Arguments

lambda

(function) the instantaneous rate of the NHPPP.

majorizer_intercept

(double) the intercept (alpha) of the loglinear majorizer function.

majorizer_slope

(double) the slope (‘beta’) of the loglinear majorizer function.

t_min

(double) the lower bound of the time interval.

t_max

(double) the upper bound of the time interval.

majorizer_is_loglinear

(boolean) if TRUE the majorizer is loglinear exp(alpha + beta * t)

atmost1

boolean, draw at most 1 event time

Value

a vector of at least 1 event times


Simulate from a zero-truncated non homogeneous Poisson Point Process (NHPPP) from (t0, t_max) (thinning method) with piecewise constant_majorizer

Description

Sample zero-truncated NHPPP times using the thinning method

Usage

ztdraw_intensity_step(lambda, majorizer_vector, time_breaks, atmost1 = FALSE)

Arguments

lambda

(function) the instantaneous rate of the NHPPP. A continuous function of time.

majorizer_vector

(scalar, double) K constant majorizing rates, one per interval

time_breaks

(vector, double) K+1 time points defining K intervals of constant rates: ⁠[t_1 = range_t[1], t_2)⁠: the first interval ⁠[t_k, t_{k+1})⁠: the k-th interval ⁠[t_{K}, t_{K+1} = range_t[2])⁠: the K-th (last) interval

atmost1

boolean, draw at most 1 event time

Value

a vector of event times (t_) with at least one element


Simulate size samples from a zero-truncated non homogeneous Poisson Point Process (zt-NHPPP) from (t_min, t_max) with linear intensity function

Description

Sample zero-truncated NHPPP times from a linear intensity function using the inversion method, optionally using an rstream generator

Usage

ztdraw_sc_linear(intercept, slope, t_min, t_max, atmost1 = FALSE)

Arguments

intercept

(double) the intercept

slope

(double) the slope

t_min

(double) the lower bound of the time interval

t_max

(double) the upper bound of the time interval

atmost1

(boolean) draw 1 event time

Value

a vector of at least 1 event times

Examples

x <- ztdraw_sc_linear(intercept = 0, slope = 0.2, t_min = 0, t_max = 10)


Simulate from a zero-truncated non homogeneous Poisson Point Process (zt-NHPPP) from (t_min, t_max) with a log-linear intensity function

Description

Sample zt-NHPPP times from an log-linear intensity function

Usage

ztdraw_sc_loglinear(intercept, slope, t_min, t_max, atmost1 = FALSE)

Arguments

intercept

(double) the intercept in the exponent

slope

(double) the slope in the exponent

t_min

(double) the lower bound of the time interval

t_max

(double) the upper bound of the time interval

atmost1

boolean, 1 event time

Value

a vector of at least 1 event times

Examples

x <- ztdraw_sc_loglinear(intercept = 0, slope = 0.2, t_min = 0, t_max = 10)


Simulate a zero-truncated homogeneous Poisson Point Process over (t_min, t_max]

Description

Simulate a zero-truncated homogeneous Poisson Point Process over (t_min, t_max]

Usage

ztppp(rate, t_min, t_max, atmost1 = FALSE)

Arguments

rate

(scalar, double) constant instantaneous rate

t_min

(scalar, double) lower bound of the time interval

t_max

(scalar, double) upper bound of the time interval

atmost1

boolean, draw at most 1 event time

Value

a vector of event times of size size

Examples

x <- ztppp(t_min = 0, t_max = 10, rate = 0.001)