Internal Documentation
Documentation for CensoredDistributions.jl's internal interface.
Contents
Index
Base.randBase.randBase.randBase.randCensoredDistributions._collect_unique_boundariesCensoredDistributions._compute_pdfs_with_cacheCensoredDistributions._gamma_cdf_ad_safeCensoredDistributions._make_weibull_gCensoredDistributions._quantile_optimizationCensoredDistributions.combine_weightsDistributions.ccdfDistributions.cdfDistributions.cdfDistributions.cdfDistributions.cdfDistributions.logccdfDistributions.logcdfDistributions.logcdfDistributions.logcdfDistributions.logcdfDistributions.logpdfDistributions.logpdfDistributions.logpdfDistributions.logpdfDistributions.logpdfDistributions.logpdfDistributions.logpdfDistributions.logpdfDistributions.pdfDistributions.pdfDistributions.pdfDistributions.pdfDistributions.pdfDistributions.samplerStatistics.meanStatistics.medianStatistics.quantileStatistics.quantileStatistics.quantileStatistics.quantileStatistics.stdStatistics.varStatsAPI.loglikelihoodStatsAPI.loglikelihoodStatsAPI.loglikelihood
Internal API
Base.rand — Methodrand(
rng::Random.AbstractRNG,
d::CensoredDistributions.IntervalCensored
) -> Any
Generate a random sample by discretising a sample from the underlying distribution.
See also: quantile
Base.rand — Methodrand(
rng::Random.AbstractRNG,
d::CensoredDistributions.PrimaryCensored
) -> Any
Generate a random sample by summing samples from delay and primary event distributions.
See also: quantile
Base.rand — Methodrand(
rng::Random.AbstractRNG,
d::CensoredDistributions.Weighted
) -> Any
Generate a random sample (delegates to underlying distribution).
See also: quantile
Base.rand — Methodrand(rng::Random.AbstractRNG, d::ExponentiallyTilted) -> Any
Generate a random sample using inverse transform sampling.
See also: quantile
CensoredDistributions._collect_unique_boundaries — Method_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.
CensoredDistributions._compute_pdfs_with_cache — Method_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.
CensoredDistributions._gamma_cdf_ad_safe — Method_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 parameterx::Real: Evaluation point
Returns
The gamma CDF value at x with shape k and scale θ.
CensoredDistributions._make_weibull_g — Method_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)CensoredDistributions._quantile_optimization — Method_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 quantilep: The probability value in[0, 1]
Keyword Arguments
initial_guess_fn: Function that takes(d, p)and returns initial guess vector. Ifnothing, usesquantile(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: Iftrue, 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)^2with 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.
CensoredDistributions.combine_weights — Methodcombine_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)
Distributions.ccdf — Methodccdf(d::CensoredDistributions.Weighted, x::Real) -> Any
Compute the complementary cumulative distribution function (delegates to underlying distribution).
See also: cdf
Distributions.cdf — Methodcdf(
d::CensoredDistributions.IntervalCensored,
x::Real
) -> Any
Compute the cumulative distribution function.
See also: logcdf
Distributions.cdf — Methodcdf(
d::CensoredDistributions.PrimaryCensored,
x::Real
) -> Any
Compute the cumulative distribution function.
See also: logcdf
Distributions.cdf — Methodcdf(d::CensoredDistributions.Weighted, x::Real) -> Any
Compute the cumulative distribution function (delegates to underlying distribution).
See also: logcdf
Distributions.cdf — Methodcdf(d::ExponentiallyTilted, x::Real) -> Any
Compute the cumulative distribution function.
Distributions.logccdf — Methodlogccdf(d::CensoredDistributions.Weighted, x::Real) -> Any
Compute the log complementary cumulative distribution function (delegates to underlying distribution).
See also: logcdf
Distributions.logcdf — Methodlogcdf(
d::CensoredDistributions.IntervalCensored,
x::Real
) -> Any
Compute the log cumulative distribution function.
See also: cdf
Distributions.logcdf — Methodlogcdf(
d::CensoredDistributions.PrimaryCensored,
x::Real
) -> Any
Compute the log cumulative distribution function.
See also: cdf
Distributions.logcdf — Methodlogcdf(d::CensoredDistributions.Weighted, x::Real) -> Any
Compute the log cumulative distribution function (delegates to underlying distribution).
See also: cdf
Distributions.logcdf — Methodlogcdf(d::ExponentiallyTilted, x::Real) -> Any
Compute the log cumulative distribution function.
See also: cdf
Distributions.logpdf — Methodlogpdf(
d::CensoredDistributions.IntervalCensored,
x::AbstractVector{<:Real}
) -> Any
Compute log probability masses for an array of values using optimised PDF computation.
Distributions.logpdf — Methodlogpdf(
d::CensoredDistributions.IntervalCensored,
x::Real
) -> Any
Compute the log probability mass for the interval containing x.
Distributions.logpdf — Methodlogpdf(
d::CensoredDistributions.PrimaryCensored,
x::Real
) -> Any
Compute the log probability density function using numerical differentiation of the log CDF.
Distributions.logpdf — Methodlogpdf(
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
Distributions.logpdf — Methodlogpdf(d::CensoredDistributions.Weighted, x::Real) -> Any
Return the weighted log-probability for scalar observations.
See also: pdf
Distributions.logpdf — Methodlogpdf(
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
Distributions.logpdf — Methodlogpdf(
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
Distributions.logpdf — Methodlogpdf(d::ExponentiallyTilted, x::Real) -> Any
Compute the log probability density function.
Distributions.pdf — Methodpdf(
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.
Distributions.pdf — Methodpdf(
d::CensoredDistributions.IntervalCensored,
x::Real
) -> Any
Compute the probability mass for the interval containing x.
Distributions.pdf — Methodpdf(
d::CensoredDistributions.PrimaryCensored,
x::Real
) -> Any
Compute the probability density function using numerical differentiation.
See also: logpdf
Distributions.pdf — Methodpdf(d::CensoredDistributions.Weighted, x::Real) -> Any
Return the probability density from the underlying distribution (unweighted).
See also: logpdf
Distributions.pdf — Methodpdf(d::ExponentiallyTilted, x::Real) -> Any
Compute the probability density function.
See also: logpdf
Distributions.sampler — Methodsampler(
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
Statistics.mean — MethodStatistics.median — Methodmedian(d::ExponentiallyTilted) -> Any
Compute the median of the distribution.
Statistics.quantile — Methodquantile(
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
Statistics.quantile — Methodquantile(
d::CensoredDistributions.PrimaryCensored,
p::Real
) -> Any
Compute the quantile (inverse CDF) using numerical optimization.
See also: cdf
Statistics.quantile — Methodquantile(d::CensoredDistributions.Weighted, p::Real) -> Any
Compute the quantile function (delegates to underlying distribution).
See also: cdf
Statistics.quantile — Methodquantile(d::ExponentiallyTilted, p::Real) -> Any
Compute the quantile function (inverse CDF) using analytical formula.
See also: cdf
Statistics.std — Methodstd(d::ExponentiallyTilted) -> Any
Compute the standard deviation of the distribution.
Statistics.var — MethodStatsAPI.loglikelihood — Methodloglikelihood(
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
StatsAPI.loglikelihood — Methodloglikelihood(
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
StatsAPI.loglikelihood — Methodloglikelihood(
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