Internal Documentation

Documentation for CensoredDistributions.jl's internal interface.

Contents

Index

Internal API

Base.randMethod
rand(
    rng::Random.AbstractRNG,
    d::CensoredDistributions.IntervalCensored
) -> Any

Generate a random sample by discretising a sample from the underlying distribution.

See also: quantile

source
Base.randMethod
rand(
    rng::Random.AbstractRNG,
    d::CensoredDistributions.PrimaryCensored
) -> Any

Generate a random sample by summing samples from delay and primary event distributions.

See also: quantile

source
Base.randMethod
rand(
    rng::Random.AbstractRNG,
    d::CensoredDistributions.Weighted
) -> Any

Generate a random sample (delegates to underlying distribution).

See also: quantile

source
Base.randMethod
rand(rng::Random.AbstractRNG, d::ExponentiallyTilted) -> Any

Generate a random sample using inverse transform sampling.

See also: quantile

source
CensoredDistributions._collect_unique_boundariesMethod
_collect_unique_boundaries(
    d::CensoredDistributions.IntervalCensored,
    x::AbstractVector{<:Real}
) -> Any
_collect_unique_boundaries(d::IntervalCensored, x::AbstractVector)

Collect all unique interval boundaries needed for vectorised PDF computation.

Returns a sorted vector of unique boundaries with appropriate type promotion.

source
CensoredDistributions._compute_pdfs_with_cacheMethod
_compute_pdfs_with_cache(
    d::CensoredDistributions.IntervalCensored,
    x::AbstractVector{<:Real},
    cdf_lookup::Dict
) -> Any
_compute_pdfs_with_cache(d::IntervalCensored, x::AbstractVector, cdf_lookup::Dict)

Compute PDFs efficiently using cached CDF values.

Uses the same boundary case handling as the scalar method.

source
CensoredDistributions._gamma_cdf_ad_safeMethod
_gamma_cdf_ad_safe(k::Real, θ::Real, x::Real) -> Any

AD-compatible gamma CDF using HypergeometricFunctions.

Based on the identity: γ(a,z) = z^a/a * M(a, a+1, -z) where γ is the lower incomplete gamma function and M is the confluent hypergeometric function.

Uses the same approach as the Weibull g function: P(a,z) = γ(a,z)/Γ(a) = z^a/a * M(a, a+1, -z) / Γ(a). For integer a, Γ(a) = (a-1)!, but uses gamma(k) for generality.

Arguments

  • k::Real: Shape parameter
  • θ::Real: Scale parameter
  • x::Real: Evaluation point

Returns

The gamma CDF value at x with shape k and scale θ.

source
CensoredDistributions._make_weibull_gMethod
_make_weibull_g(
    k::Real,
    λ::Real
) -> CensoredDistributions.var"#weibull_g_specialized#1"{<:Real, <:Real}

Function factory for optimized Weibull g function.

Creates a specialized function with pre-computed constants for g(t; k, λ) = γ(1 + 1/k, (t/λ)^k) where γ is the lower incomplete gamma function.

Pre-computes inv_k = 1/k and a = 1 + inv_k to avoid repeated computation in the returned specialized function.

Arguments

  • k::Real: Weibull shape parameter
  • λ::Real: Weibull scale parameter

Returns

A specialized function weibull_g_specialized(t::Real) that efficiently computes the Weibull g function using pre-computed constants.

Examples

weibull_g_func = _make_weibull_g(2.0, 1.5)
g_val = weibull_g_func(3.0)
source
CensoredDistributions._quantile_optimizationMethod
_quantile_optimization(
    d,
    p::Real;
    initial_guess_fn,
    result_postprocess_fn,
    check_nan
) -> Any

Internal function for quantile optimization using numerical methods.

Solves the equation cdf(d, q) - p = 0 using the Nelder-Mead algorithm. This is shared logic used by both PrimaryCensored and IntervalCensored quantile functions.

Arguments

  • d: The distribution for which to compute the quantile
  • p: The probability value in [0, 1]

Keyword Arguments

  • initial_guess_fn: Function that takes (d, p) and returns initial guess vector. If nothing, uses quantile(get_dist(d), p) as scalar initial guess.
  • result_postprocess_fn: Function to post-process the optimization result. Defaults to identity (no post-processing).
  • check_nan: If true, explicitly check for NaN input values.

Returns

The quantile value after optimization and post-processing.

Implementation Details

  • Validates that p ∈ [0, 1] (with optional NaN checking)
  • Handles boundary cases p=0 (minimum) and p=1 (maximum) analytically
  • Creates objective function (cdf(d, q) - p)^2 with support checking
  • Uses heavy penalty for values outside distribution support
  • Solves with Nelder-Mead algorithm with tight tolerances
  • Checks convergence and applies post-processing to result

Type Stability

This function is designed to maintain type stability when the initial guess function and post-processing function are type-stable.

source
CensoredDistributions.combine_weightsMethod
combine_weights(_::Missing, _::Missing) -> Missing

Combine constructor weight with observation weight using dispatch-based rules.

Weight combination rules:

  • missing, missing → missing (both missing means no weight)
  • w1, missing → w1 (use constructor weight)
  • missing, w2 → w2 (use observation weight)
  • w1, w2 → w1 * w2 (multiply weights)

Vector Extensions

For Product distributions, additional methods handle vectorised weight combinations:

  • Vector, Vector → combine_weights.(vector1, vector2) (element-wise combination)
  • Vector, missing → Vector (keep constructor weights)
  • Vector, scalar → [combine_weights(w, scalar) for w in Vector] (broadcast scalar)
source
Distributions.ccdfMethod
ccdf(d::CensoredDistributions.Weighted, x::Real) -> Any

Compute the complementary cumulative distribution function (delegates to underlying distribution).

See also: cdf

source
Distributions.cdfMethod
cdf(
    d::CensoredDistributions.IntervalCensored,
    x::Real
) -> Any

Compute the cumulative distribution function.

See also: logcdf

source
Distributions.cdfMethod
cdf(
    d::CensoredDistributions.PrimaryCensored,
    x::Real
) -> Any

Compute the cumulative distribution function.

See also: logcdf

source
Distributions.cdfMethod
cdf(d::CensoredDistributions.Weighted, x::Real) -> Any

Compute the cumulative distribution function (delegates to underlying distribution).

See also: logcdf

source
Distributions.logccdfMethod
logccdf(d::CensoredDistributions.Weighted, x::Real) -> Any

Compute the log complementary cumulative distribution function (delegates to underlying distribution).

See also: logcdf

source
Distributions.logcdfMethod
logcdf(
    d::CensoredDistributions.IntervalCensored,
    x::Real
) -> Any

Compute the log cumulative distribution function.

See also: cdf

source
Distributions.logcdfMethod
logcdf(
    d::CensoredDistributions.PrimaryCensored,
    x::Real
) -> Any

Compute the log cumulative distribution function.

See also: cdf

source
Distributions.logcdfMethod
logcdf(d::CensoredDistributions.Weighted, x::Real) -> Any

Compute the log cumulative distribution function (delegates to underlying distribution).

See also: cdf

source
Distributions.logcdfMethod
logcdf(d::ExponentiallyTilted, x::Real) -> Any

Compute the log cumulative distribution function.

See also: cdf

source
Distributions.logpdfMethod
logpdf(
    d::CensoredDistributions.IntervalCensored,
    x::AbstractVector{<:Real}
) -> Any

Compute log probability masses for an array of values using optimised PDF computation.

See also: pdf, logpdf

source
Distributions.logpdfMethod
logpdf(
    d::CensoredDistributions.IntervalCensored,
    x::Real
) -> Any

Compute the log probability mass for the interval containing x.

See also: pdf, logcdf

source
Distributions.logpdfMethod
logpdf(
    d::CensoredDistributions.PrimaryCensored,
    x::Real
) -> Any

Compute the log probability density function using numerical differentiation of the log CDF.

See also: pdf, logcdf

source
Distributions.logpdfMethod
logpdf(
    d::CensoredDistributions.Weighted,
    obs::NamedTuple{(:value, :weight)}
) -> Any

Return the weighted log-probability for joint observations as NamedTuple.

Combines constructor weight with observation weight via multiplication. Expected format: (value = x, weight = w).

See also: pdf

source
Distributions.logpdfMethod
logpdf(d::CensoredDistributions.Weighted, x::Real) -> Any

Return the weighted log-probability for scalar observations.

See also: pdf

source
Distributions.logpdfMethod
logpdf(
    d::Distributions.Product{<:Distributions.ValueSupport, <:CensoredDistributions.Weighted, <:AbstractVector{<:CensoredDistributions.Weighted}},
    x::AbstractVector{<:Real}
) -> Any

Efficient vectorised log-probability computation for Product{<:ValueSupport, <:Weighted} with vector observations.

See also: logpdf

source
Distributions.logpdfMethod
logpdf(
    d::Distributions.Product{<:Distributions.ValueSupport, <:CensoredDistributions.Weighted, <:AbstractVector{<:CensoredDistributions.Weighted}},
    obs::NamedTuple{(:values, :weights)}
) -> Any

Efficient vectorised log-probability computation for Product{<:ValueSupport, <:Weighted} with joint observations.

Handles joint observations and weight stacking. Expected format: (values = [...], weights = [...]).

See also: logpdf

source
Distributions.pdfMethod
pdf(
    d::CensoredDistributions.IntervalCensored,
    x::AbstractVector{<:Real}
) -> Any

Compute probability masses for an array of values using optimised vectorisation.

This method collects unique interval boundaries, computes CDFs once, then uses cached values for efficient PDF computation across the array.

See also: pdf, logpdf

source
Distributions.pdfMethod
pdf(
    d::CensoredDistributions.IntervalCensored,
    x::Real
) -> Any

Compute the probability mass for the interval containing x.

See also: logpdf, cdf

source
Distributions.pdfMethod
pdf(
    d::CensoredDistributions.PrimaryCensored,
    x::Real
) -> Any

Compute the probability density function using numerical differentiation.

See also: logpdf

source
Distributions.pdfMethod
pdf(d::CensoredDistributions.Weighted, x::Real) -> Any

Return the probability density from the underlying distribution (unweighted).

See also: logpdf

source
Distributions.samplerMethod
sampler(
    d::CensoredDistributions.Weighted
) -> Union{CensoredDistributions.Weighted{D, Missing} where D<:(Distributions.UnivariateDistribution), CensoredDistributions.Weighted{D, T} where {D<:(Distributions.UnivariateDistribution), T<:Real}}

Create a sampler for efficient sampling (delegates to underlying distribution).

See also: rand

source
Statistics.quantileMethod
quantile(
    d::CensoredDistributions.IntervalCensored,
    p::Real
) -> Any

Compute the quantile using numerical optimization.

The returned quantile respects the interval structure:

  • For regular intervals: quantiles are multiples of the interval width
  • For arbitrary intervals: quantiles correspond to interval boundary values

See also: cdf

source
Statistics.quantileMethod
quantile(
    d::CensoredDistributions.PrimaryCensored,
    p::Real
) -> Any

Compute the quantile (inverse CDF) using numerical optimization.

See also: cdf

source
Statistics.quantileMethod
quantile(d::CensoredDistributions.Weighted, p::Real) -> Any

Compute the quantile function (delegates to underlying distribution).

See also: cdf

source
Statistics.quantileMethod
quantile(d::ExponentiallyTilted, p::Real) -> Any

Compute the quantile function (inverse CDF) using analytical formula.

See also: cdf

source
Statistics.stdMethod
std(d::ExponentiallyTilted) -> Any

Compute the standard deviation of the distribution.

See also: mean, var

source
StatsAPI.loglikelihoodMethod
loglikelihood(
    d::CensoredDistributions.Weighted,
    obs::NamedTuple{(:value, :weight)}
) -> Any

Compute log-likelihood for single Weighted distribution with joint observations.

Handles joint observations as NamedTuple: (value = x, weight = w).

See also: logpdf

source
StatsAPI.loglikelihoodMethod
loglikelihood(
    d::CensoredDistributions.Weighted,
    obs::NamedTuple{(:values, :weights)}
) -> Any

Compute log-likelihood for single Weighted distribution with vectorized joint observations.

Handles joint observations as NamedTuple: (values = [...], weights = [...]). This is useful when a single weighted distribution is used with multiple observations.

See also: logpdf

source
StatsAPI.loglikelihoodMethod
loglikelihood(
    d::Distributions.Product{<:Distributions.ValueSupport, <:CensoredDistributions.Weighted, <:AbstractVector{<:CensoredDistributions.Weighted}},
    obs::NamedTuple{(:values, :weights)}
) -> Any

Compute log-likelihood for Product{<:ValueSupport, <:Weighted} with joint observations.

Handles joint observations as NamedTuple: (values = [...], weights = [...]).

See also: logpdf

source