Public Documentation
Documentation for CensoredDistributions.jl
's public interface.
See the Internals section of the manual for internal package docs covering all submodules.
Contents
Index
CensoredDistributions.AbstractSolverMethod
CensoredDistributions.AnalyticalSolver
CensoredDistributions.ExponentiallyTilted
CensoredDistributions.ExponentiallyTilted
CensoredDistributions.IntervalCensored
CensoredDistributions.NumericSolver
CensoredDistributions.PrimaryCensored
CensoredDistributions.Weighted
CensoredDistributions.double_interval_censored
CensoredDistributions.get_dist
CensoredDistributions.get_dist
CensoredDistributions.get_dist
CensoredDistributions.get_dist
CensoredDistributions.get_dist
CensoredDistributions.get_dist
CensoredDistributions.get_dist
CensoredDistributions.get_dist_recursive
CensoredDistributions.interval_censored
CensoredDistributions.interval_censored
CensoredDistributions.primary_censored
CensoredDistributions.primary_censored
CensoredDistributions.primarycensored_cdf
CensoredDistributions.primarycensored_cdf
CensoredDistributions.primarycensored_cdf
CensoredDistributions.primarycensored_cdf
CensoredDistributions.primarycensored_cdf
CensoredDistributions.primarycensored_cdf
CensoredDistributions.primarycensored_logcdf
CensoredDistributions.weight
CensoredDistributions.weight
CensoredDistributions.weight
CensoredDistributions.weight
CensoredDistributions.weight
Public API
CensoredDistributions.AbstractSolverMethod
— Typeabstract type AbstractSolverMethod
Abstract type for solver methods used in CDF computation.
Subtypes determine whether analytical solutions are preferred or numerical integration is forced.
Fields
CensoredDistributions.AnalyticalSolver
— Typestruct AnalyticalSolver{S} <: CensoredDistributions.AbstractSolverMethod
Solver that attempts analytical solutions when available, falling back to numerical integration.
Stores a numerical integration solver for use when no analytical solution exists for a given distribution pair.
Fields
solver::Any
CensoredDistributions.ExponentiallyTilted
— Typestruct ExponentiallyTilted{T<:Real} <: Distributions.Distribution{Distributions.Univariate, Distributions.Continuous}
A continuous distribution on interval [min, max] with exponential tilting controlled by parameter r. This distribution generalises the uniform distribution by allowing exponential weighting of values within the interval.
Mathematical Definition
The probability density function is:
\[f(x) = \frac{r \exp(r(x - \text{min}))}{(\exp(r(\text{max} - \text{min})) - 1)}\]
for x ∈ [min, max] and r ≠ 0.
The cumulative distribution function is:
\[F(x) = \frac{\exp(r(x - \text{min})) - 1}{\exp(r(\text{max} - \text{min})) - 1}\]
for x ∈ [min, max].
The quantile function (inverse CDF) is:
\[F^{-1}(p) = \text{min} + \frac{1}{r} \log\left(1 + p(\exp(r(\text{max} - \text{min})) - 1)\right)\]
for p ∈ [0, 1].
When r → 0, all functions reduce to the uniform distribution on [min, max].
- For r > 0: distribution is tilted towards higher values (increasing density)
- For r < 0: distribution is tilted towards lower values (decreasing density)
Fields
min::Real
: Lower bound of the distribution support.max::Real
: Upper bound of the distribution support.r::Real
: Growth rate parameter (tilting parameter).
CensoredDistributions.ExponentiallyTilted
— MethodCreate an exponentially tilted distribution.
An exponentially tilted distribution on interval [min, max] with growth rate r. For r > 0, the distribution is tilted towards higher values; for r < 0, towards lower values. When r ≈ 0, it approaches a uniform distribution.
Examples
using CensoredDistributions, Distributions
# Exponentially increasing distribution
d1 = ExponentiallyTilted(0.0, 1.0, 2.0)
pdf(d1, 0.5)
# Nearly uniform distribution (small r)
d3 = ExponentiallyTilted(0.0, 1.0, 1e-10)
rand(d3, 5)
CensoredDistributions.IntervalCensored
— Typestruct IntervalCensored{D<:(Distributions.UnivariateDistribution), T} <: Distributions.UnivariateDistribution{Distributions.ValueSupport}
Interval-censored distribution where continuous values are observed only within intervals.
Supports both:
- Regular intervals: fixed-width intervals (e.g., daily reporting)
- Arbitrary intervals: custom interval boundaries
Fields
dist::Distributions.UnivariateDistribution
: The underlying continuous distributionboundaries::Any
: Either a scalar (regular intervals) or vector (boundaries for arbitrary intervals)
CensoredDistributions.NumericSolver
— Typestruct NumericSolver{S} <: CensoredDistributions.AbstractSolverMethod
Solver that always uses numerical integration.
Forces numerical computation even when analytical solutions are available, useful for testing and validation.
The solver
field contains the numerical integration solver to use.
Fields
solver::Any
CensoredDistributions.PrimaryCensored
— Typestruct PrimaryCensored{D1<:(Distributions.UnivariateDistribution), D2<:(Distributions.UnivariateDistribution), M<:CensoredDistributions.AbstractSolverMethod} <: Distributions.Distribution{Distributions.Univariate, Distributions.Continuous}
Represents the distribution of observed delays when the primary event time is subject to censoring.
The dist
field contains the delay distribution from primary event to observation. The primary_event
field contains the primary event time distribution. The method
field determines computation strategy:
AnalyticalSolver
: Uses closed-form solutions when available (Gamma, LogNormal, Weibull with Uniform primary), falls back to numeric otherwiseNumericSolver
: Always uses quadrature integration
See also
primary_censored
: Constructor functionprimarycensored_cdf
: CDF computation with method dispatch
Fields
dist::Distributions.UnivariateDistribution
: The delay distribution from primary event to observation.primary_event::Distributions.UnivariateDistribution
: The primary event time distribution.method::CensoredDistributions.AbstractSolverMethod
: The solver method for CDF computation.
CensoredDistributions.Weighted
— Typestruct Weighted{D<:(Distributions.UnivariateDistribution), T<:Union{Missing, Real}} <: Distributions.UnivariateDistribution{Distributions.ValueSupport}
A distribution wrapper that applies a weight to the log-probability of an underlying distribution. This is primarily used where observations have associated counts or weights.
Only the logpdf
method is affected by the weight - all other methods (pdf, cdf, sampling, etc.) delegate directly to the underlying distribution.
Weight Types Supported
The Weighted
struct supports three different weight scenarios:
- Real weights: Constructor weight is a specific value (e.g.,
2.5
) - Missing weights: Constructor weight is
missing
, allowing weights to be provided at observation time via joint observations(value = x, weight = w)
- Zero weights: Handled specially to return
-Inf
and avoid NaN from0 * -Inf
Examples
using CensoredDistributions, Distributions
# Single weighted observation
d = LogNormal(1.5, 0.5)
wd = weight(d, 10.0) # Observation with weight/count of 10
# Weighted log-probability calculation
observed_value = 2.0
weighted_logpdf = logpdf(wd, observed_value)
# Compare with manual calculation
manual_logpdf = 10.0 * logpdf(d, observed_value)
# weighted_logpdf ≈ manual_logpdf
Fields
dist::Distributions.UnivariateDistribution
: The underlying distribution being weighted.weight::Union{Missing, Real}
: The weight to apply to log-probabilities.
CensoredDistributions.double_interval_censored
— Methoddouble_interval_censored(
dist::Distributions.UnivariateDistribution;
primary_event,
lower,
upper,
interval,
force_numeric
) -> Union{CensoredDistributions.PrimaryCensored{D1, Distributions.Uniform{Float64}, CensoredDistributions.AnalyticalSolver{Integrals.QuadGKJL{typeof(LinearAlgebra.norm), Nothing}}} where D1<:(Distributions.UnivariateDistribution), CensoredDistributions.PrimaryCensored{D1, Distributions.Uniform{Float64}, CensoredDistributions.NumericSolver{Integrals.QuadGKJL{typeof(LinearAlgebra.norm), Nothing}}} where D1<:(Distributions.UnivariateDistribution)}
Create a distribution that combines primary interval censoring, optional truncation, and optional secondary interval censoring in the correct order.
This is a convenience function that applies the following transformations in sequence:
- Primary event censoring using
primary_censored(dist, primary_event)
- Truncation using
truncated(dist; lower=lower, upper=upper)
(if bounds are provided) - Secondary interval censoring using
interval_censored(dist, interval)
(if interval is provided)
The order of operations ensures mathematical correctness, particularly that truncation occurs before secondary interval censoring.
Arguments
dist
: The delay distribution from primary event to observation
Keyword Arguments
primary_event
: The primary event time distribution. Defaults to uniform distribution over [0, 1].lower
: Lower truncation bound. Ifnothing
, no lower truncation is applied.upper
: Upper truncation bound (e.g., observation timeD
). Ifnothing
, no upper truncation is applied.interval
: Secondary censoring interval width (e.g., daily reporting). Ifnothing
, no interval censoring is applied.force_numeric
: Force numerical integration for primary censoring even when analytical solutions exist.
Returns
A composed distribution that can be used with all standard Distributions.jl
methods (rand
, pdf
, cdf
, etc.).
Examples
using CensoredDistributions, Distributions
# Basic primary censoring only (uses default Uniform(0, 1))
dist1 = double_interval_censored(LogNormal(1.5, 0.75))
# Primary censoring + truncation
dist2 = double_interval_censored(LogNormal(1.5, 0.75); upper=10)
# Primary censoring + secondary interval censoring
dist3 = double_interval_censored(LogNormal(1.5, 0.75); interval=1)
# Full double interval censoring with truncation
dist4 = double_interval_censored(LogNormal(1.5, 0.75); upper=10, interval=1)
# Evaluate distribution functions and compute statistics
pdf_at_3 = pdf(dist4, 3.0) # probability density at 3 days
cdf_at_7 = cdf(dist4, 7.0) # P(X ≤ 7) with full censoring pipeline
ccdf_at_5 = ccdf(dist4, 5.0) # survival function
q25 = quantile(dist4, 0.25) # 25th percentile
median_delay = quantile(dist4, 0.5) # median
q95 = quantile(dist4, 0.95) # 95th percentile
samples = rand(dist4, 100) # random samples
# Custom primary event distribution
dist5 = double_interval_censored(LogNormal(1.5, 0.75); primary_event=Uniform(0, 2))
Mathematical Background
This function implements the complete workflow for handling censored delay distributions as described in Park et al. (2024) and Charniga et al. (2024):
- Primary censoring: Accounts for uncertainty in the primary event timing
- Truncation: Handles observation windows and finite study periods
- Secondary censoring: Models interval censoring effects (e.g., daily reporting)
References
- Park et al. (2024): "Estimating epidemiological delay distributions for infectious diseases"
- Charniga et al. (2024): "Best practices for estimating and reporting epidemiological delay distributions"
CensoredDistributions.get_dist
— Methodget_dist(d) -> Distributions.UnivariateDistribution
Extract the underlying distribution from a wrapped distribution type.
This utility function provides a consistent interface for extracting the core distribution from various wrapper types in the CensoredDistributions ecosystem. For unwrapped distributions, it returns the distribution unchanged.
Arguments
d
: A distribution or wrapped distribution
Returns
The underlying distribution. For base distributions, returns d
unchanged.
Examples
using CensoredDistributions, Distributions
# Base distribution - returns unchanged
d1 = Normal(0, 1)
get_dist(d1) == d1 # true
# Primary censored distribution
delay = LogNormal(1.5, 0.75)
window = Uniform(0, 1)
pc = primary_censored(delay, window)
get_dist(pc) == delay # true
# Interval censored distribution
continuous = Normal(5, 2)
ic = interval_censored(continuous, 1.0)
get_dist(ic) == continuous # true
CensoredDistributions.get_dist
— Methodget_dist(
d::CensoredDistributions.IntervalCensored
) -> Distributions.UnivariateDistribution
Extract the underlying continuous distribution from an interval censored distribution.
Returns the continuous distribution that was discretised into intervals.
CensoredDistributions.get_dist
— Methodget_dist(
d::CensoredDistributions.PrimaryCensored
) -> Distributions.UnivariateDistribution
Extract the delay distribution from a primary censored distribution.
Returns the underlying delay distribution (the distribution of times from primary event to observation).
CensoredDistributions.get_dist
— Methodget_dist(
d::CensoredDistributions.Weighted
) -> Distributions.UnivariateDistribution
Extract the underlying distribution from a weighted distribution.
Returns the base distribution before weighting was applied.
CensoredDistributions.get_dist
— Methodget_dist(
d::Distributions.Censored
) -> Distributions.UnivariateDistribution
Extract the underlying distribution from a censored distribution.
Returns the uncensored distribution before censoring bounds were applied. This method supports the Censored
type from Distributions.jl created by the censored()
function.
Examples
using Distributions
# Extract Normal from censored Normal
base = Normal(0, 1)
censored_dist = censored(base, -2, 2)
get_dist(censored_dist) === base # true
# Works with any censored distribution
gamma_base = Gamma(2, 1)
censored_gamma = censored(gamma_base, 0.1, 5.0)
get_dist(censored_gamma) === gamma_base # true
CensoredDistributions.get_dist
— Methodget_dist(
d::Distributions.Product
) -> AbstractVector{T} where {S<:Distributions.ValueSupport, T<:Distributions.UnivariateDistribution{S}}
Extract the component distributions from a product distribution.
Returns a vector containing the underlying distributions for each component of the product distribution. This is useful when working with vectorised distributions or multiple independent observations.
Note
This method has different behaviour from other get_dist
methods as it returns a vector of distributions rather than a single distribution.
CensoredDistributions.get_dist
— Methodget_dist(
d::Distributions.Truncated
) -> Distributions.UnivariateDistribution
Extract the untruncated distribution from a truncated distribution.
Returns the underlying continuous distribution before truncation bounds were applied.
Examples
using Distributions
# Extract Normal from truncated Normal
base = Normal(0, 1)
trunc_dist = truncated(base, -2, 2)
get_dist(trunc_dist) === base # true
# Works with any truncated distribution
gamma_base = Gamma(2, 1)
trunc_gamma = truncated(gamma_base, 0.1, 5.0)
get_dist(trunc_gamma) === gamma_base # true
CensoredDistributions.get_dist_recursive
— Methodget_dist_recursive(d) -> Any
Recursively extract the underlying distribution from nested wrapper types.
This function keeps applying get_dist
until it reaches a distribution that doesn't have a specialised method, meaning no further unwrapping is possible. This is useful for deeply nested distributions like IntervalCensored{Truncated{PrimaryCensored{...}}}
.
Arguments
d
: A distribution or nested wrapped distribution
Returns
The deeply underlying distribution after all unwrapping is complete.
Examples
using CensoredDistributions, Distributions
# Single wrapper - same as get_dist
delay = LogNormal(1.5, 0.75)
window = Uniform(0, 1)
pc = primary_censored(delay, window)
get_dist_recursive(pc) == delay # true
# Nested wrappers
continuous = Normal(5, 2)
ic = interval_censored(continuous, 1.0)
weighted = weight(ic, 2.0)
get_dist_recursive(weighted) == continuous # true
# Truncated nested wrappers
base = Normal(0, 1)
trunc_dist = truncated(base, -2, 2)
weighted_trunc = weight(trunc_dist, 2.0)
get_dist_recursive(weighted_trunc) == base # true
# Censored nested wrappers
censored_dist = censored(Normal(0, 1), -2, 2)
weighted_censored = weight(censored_dist, 2.0)
get_dist_recursive(weighted_censored) == base # true
# Double interval censored distributions (fully recursive)
delay = LogNormal(1.5, 0.75)
full_double = double_interval_censored(delay; upper=10, interval=1)
get_dist_recursive(full_double) == delay # true
# Weighted double interval censored
weighted_double = weight(full_double, 2.0)
get_dist_recursive(weighted_double) == delay # true
# Base distribution - returns unchanged
d = Normal(0, 1)
get_dist_recursive(d) == d # true
Note
For Product
distributions, this function applies recursive extraction to each component, potentially returning mixed types of underlying distributions.
CensoredDistributions.interval_censored
— Methodinterval_censored(
dist::Distributions.UnivariateDistribution,
boundaries::AbstractVector{<:Real}
) -> CensoredDistributions.IntervalCensored{D, <:AbstractVector{var"#s11"}} where {D<:(Distributions.UnivariateDistribution), var"#s11"<:Real}
Construct an interval-censored distribution with arbitrary intervals.
Creates a distribution where observations are censored to specified intervals defined by the boundaries. For example, with boundaries=[0, 2, 5, 10]
, observations fall into [0,2), [2,5), or [5,10).
Arguments
dist
: The underlying continuous distributionboundaries
: Vector of interval boundaries (must be sorted and strictly increasing, minimum length 2)
Examples
using CensoredDistributions, Distributions
# Age groups: 0-18, 18-65, 65+
age_dist = interval_censored(Normal(40, 20), [0, 18, 65, 100])
q75 = quantile(age_dist, 0.75) # 75th percentile (boundary value)
CensoredDistributions.interval_censored
— Methodinterval_censored(
dist::Distributions.UnivariateDistribution,
interval::Real
) -> CensoredDistributions.IntervalCensored{D, <:Real} where D<:(Distributions.UnivariateDistribution)
Construct an interval-censored distribution with regular intervals.
Creates a distribution where observations are censored to regular intervals of width interval
, starting from 0. For example, with interval=1
, observations fall into [0,1), [1,2), [2,3), etc.
Arguments
dist
: The underlying continuous distributioninterval
: Width of regular intervals (must be positive)
Examples
using CensoredDistributions, Distributions
# Daily reporting intervals
d = interval_censored(Normal(5, 2), 1.0)
# Evaluate distribution functions
pdf_at_5 = pdf(d, 5.0) # probability mass at interval containing 5
q25 = quantile(d, 0.25) # 25th percentile (interval boundary)
CensoredDistributions.primary_censored
— Methodprimary_censored(
dist::Distributions.UnivariateDistribution,
primary_event::Distributions.UnivariateDistribution;
solver,
force_numeric
) -> Union{CensoredDistributions.PrimaryCensored{D1, D2, CensoredDistributions.AnalyticalSolver{Integrals.QuadGKJL{typeof(LinearAlgebra.norm), Nothing}}} where {D1<:(Distributions.UnivariateDistribution), D2<:(Distributions.UnivariateDistribution)}, CensoredDistributions.PrimaryCensored{D1, D2, CensoredDistributions.NumericSolver{Integrals.QuadGKJL{typeof(LinearAlgebra.norm), Nothing}}} where {D1<:(Distributions.UnivariateDistribution), D2<:(Distributions.UnivariateDistribution)}}
Create a primary event censored distribution.
Models a process where a primary event occurs within a censoring window, followed by a delay. The primary event time is not observed directly but is known to fall within the censoring distribution's support. The observed time is the sum of the primary event time and the delay.
Method Selection
The CDF computation is handled by primarycensored_cdf
which automatically dispatches between:
- Analytical methods: Available for these distribution pairs with Uniform primary events when
force_numeric=false
:Gamma
delay distributionLogNormal
delay distributionWeibull
delay distribution
- Numeric integration: Falls back to quadrature for all other distribution pairs or when
force_numeric=true
Set force_numeric=true
to always use numeric integration, which may be necessary for certain AD backends or when debugging.
Arguments
dist
: The delay distribution from primary event to observationprimary_event
: The distribution of primary event times within the window
Keyword Arguments
solver
: Quadrature solver (default:QuadGKJL()
)force_numeric
: Force numeric integration even when analytical available (default:false
)
This is useful for modeling:
- Infection-to-symptom onset times when infection time is uncertain
- Exposure-to-outcome delays with uncertain exposure timing
- Any process where the initiating event time has uncertainty
Examples
using CensoredDistributions, Distributions
# Incubation period (delay) with uncertain infection time (primary event)
incubation = LogNormal(1.5, 0.75) # Delay distribution
infection_window = Uniform(0, 1) # Daily infection window
d = primary_censored(incubation, infection_window)
# Evaluate distribution functions
pdf_at_2 = pdf(d, 2.0) # probability density at 2 days
cdf_at_5 = cdf(d, 5.0) # cumulative probability by 5 days
q50 = quantile(d, 0.5) # median
# Force numeric integration for debugging or AD compatibility
d_numeric = primary_censored(incubation, infection_window; force_numeric=true)
See also
primarycensored_cdf
: The underlying CDF computation with method dispatch
CensoredDistributions.primary_censored
— Methodprimary_censored(
dist::Distributions.UnivariateDistribution;
primary_event,
solver,
force_numeric
) -> Union{CensoredDistributions.PrimaryCensored{D1, Distributions.Uniform{Float64}, CensoredDistributions.AnalyticalSolver{Integrals.QuadGKJL{typeof(LinearAlgebra.norm), Nothing}}} where D1<:(Distributions.UnivariateDistribution), CensoredDistributions.PrimaryCensored{D1, Distributions.Uniform{Float64}, CensoredDistributions.NumericSolver{Integrals.QuadGKJL{typeof(LinearAlgebra.norm), Nothing}}} where D1<:(Distributions.UnivariateDistribution)}
Create a primary event censored distribution with keyword arguments.
This is a convenience version of primary_censored
that uses keyword arguments for consistency with double_interval_censored
. The primary event distribution defaults to Uniform(0, 1)
.
Examples
using CensoredDistributions, Distributions
# Using default Uniform(0, 1) primary event
d1 = primary_censored(LogNormal(1.5, 0.75))
# Custom primary event distribution
d2 = primary_censored(LogNormal(1.5, 0.75); primary_event=Uniform(0, 2))
CensoredDistributions.primarycensored_cdf
— Methodprimarycensored_cdf(
dist::Distributions.Gamma,
primary_event::Distributions.Uniform,
x::Real,
_::CensoredDistributions.AnalyticalSolver
) -> Any
Analytical CDF for Gamma delay with Uniform primary event distribution.
Implements the closed-form solution derived in Park et al. (2024) and originally in Cori et al. (2013). Uses the partial expectation of the Gamma distribution to avoid numerical integration.
The formula involves CDFs of Gamma distributions with shape parameters k and k+1, computed in log-space for numerical stability.
CensoredDistributions.primarycensored_cdf
— Methodprimarycensored_cdf(
dist::Distributions.LogNormal,
primary_event::Distributions.Uniform,
x::Real,
_::CensoredDistributions.AnalyticalSolver
) -> Any
Analytical CDF for LogNormal delay with Uniform primary event distribution.
Uses a parameter shift approach where the partial expectation of LogNormal(μ, σ) can be expressed using the CDF of LogNormal(μ + σ², σ).
CensoredDistributions.primarycensored_cdf
— Methodprimarycensored_cdf(
dist::Distributions.Weibull,
primary_event::Distributions.Uniform,
x::Real,
_::CensoredDistributions.AnalyticalSolver
) -> Any
Analytical CDF for Weibull delay with Uniform primary event distribution.
Uses the lower incomplete gamma function to express the partial expectation of the Weibull distribution analytically.
CensoredDistributions.primarycensored_cdf
— Methodprimarycensored_cdf(
dist::Distributions.UnivariateDistribution,
primary_event::Distributions.UnivariateDistribution,
x::Real,
method::CensoredDistributions.AbstractSolverMethod
) -> Any
Compute the CDF of a primary event censored distribution.
Dispatches to either analytical or numerical implementation based on the solver method. Analytical solutions are available for specific distribution pairs with Uniform primary events. Use methods(primarycensored_cdf)
to see all available analytical implementations.
For distribution pairs without analytical solutions, numerical integration is used automatically when called with an AnalyticalSolver
. Use NumericSolver
to force numerical integration even when analytical solutions exist (useful for testing and validation).
Arguments
dist
: The delay distribution from primary event to observationprimary_event
: The primary event time distributionx
: Evaluation point for the CDFmethod
: Solver method (AnalyticalSolver
orNumericSolver
)
Returns
The cumulative probability P(X ≤ x) where X is the observed delay time.
Examples
using CensoredDistributions, Distributions, Integrals
# Analytical solution for Gamma with Uniform primary event
gamma_dist = Gamma(2.0, 1.5)
primary_uniform = Uniform(0, 2)
analytical_method = AnalyticalSolver(QuadGKJL())
cdf_val = primarycensored_cdf(gamma_dist, primary_uniform, 3.0, analytical_method)
# Force numerical integration
numeric_method = NumericSolver(QuadGKJL())
cdf_numeric = primarycensored_cdf(gamma_dist, primary_uniform, 3.0, numeric_method)
CensoredDistributions.primarycensored_cdf
— Methodprimarycensored_cdf(
dist::Distributions.UnivariateDistribution,
primary_event::Distributions.UnivariateDistribution,
x::Real,
method::CensoredDistributions.AnalyticalSolver
) -> Any
Generic fallback implementation for AnalyticalSolver.
When no specific analytical solution is available for a distribution pair, this method falls back to numerical integration using the solver stored in the AnalyticalSolver.
Specific analytical implementations exist for:
- Gamma delay with Uniform primary event
- LogNormal delay with Uniform primary event
- Weibull delay with Uniform primary event
For all other distribution pairs, numerical integration is used automatically.
CensoredDistributions.primarycensored_cdf
— Methodprimarycensored_cdf(
dist::Distributions.UnivariateDistribution,
primary_event::Distributions.UnivariateDistribution,
x::Real,
method::CensoredDistributions.NumericSolver
) -> Any
Numerical CDF implementation for primary event censored distributions.
Computes the CDF using numerical integration when no analytical solution is available. The integral computed is:
\[F_{S+}(x) = \int_{\max(x-u_{\max}, s_{\min})}^{x-u_{\min}} F_S(u) f_U(x-u) du\]
where FS is the delay distribution CDF, fU is the primary event distribution PDF, umin and umax are the primary event bounds, and s_min is the delay minimum.
Handles edge cases and uses the solver stored in the NumericSolver method.
CensoredDistributions.primarycensored_logcdf
— Methodprimarycensored_logcdf(
dist::Distributions.UnivariateDistribution,
primary_event::Distributions.UnivariateDistribution,
x::Real,
method::CensoredDistributions.AbstractSolverMethod
) -> Any
Compute the log CDF of a primary event censored distribution.
Dispatches to either analytical or numerical implementation based on the solver method. Computes log(CDF) with appropriate handling for edge cases where CDF = 0. Returns -Inf when the probability is zero or when numerical issues occur.
Arguments
dist
: The delay distribution from primary event to observationprimary_event
: The primary event time distributionx
: Evaluation point for the log CDFmethod
: Solver method (AnalyticalSolver
orNumericSolver
)
Examples
using CensoredDistributions, Distributions
# Create a Gamma delay with uniform primary event
dist = Gamma(2.0, 1.5)
primary_event = Uniform(0.0, 3.0)
method = AnalyticalSolver(nothing)
# Compute log CDF at x = 5.0
log_prob = primarycensored_logcdf(dist, primary_event, 5.0, method)
See also
primarycensored_cdf
: The underlying CDF computationcdf
: Interface method for PrimaryCensored distributionslogcdf
: Interface method for PrimaryCensored distributions
CensoredDistributions.weight
— Methodweight(
dists::AbstractVector{<:Distributions.UnivariateDistribution},
weights::AbstractVector{<:Real}
) -> Any
Create a product distribution of weighted distributions, where each distribution has its own weight.
A Product
distribution of Weighted
distributions suitable for vectorized observations with different distributions.
Examples
using CensoredDistributions, Distributions
y_obs = [3.5, 4.2, 3.8] # Observed values
dists = [Normal(2.0, 0.5), Normal(2.5, 0.8), Normal(1.8, 0.6)]
n_counts = [25, 10, 15] # Counts for each observation
weighted_dists = weight(dists, n_counts)
CensoredDistributions.weight
— Methodweight(
dists::AbstractVector{<:Distributions.UnivariateDistribution}
) -> Any
Create a product distribution of weighted distributions with missing constructor weights.
Useful for creating distributions where weights will be provided at observation time. Each distribution uses missing
as constructor weight, enabling observation weight to be used directly.
Examples
using CensoredDistributions, Distributions
y_obs = [3.5, 4.2, 3.8] # Observed values
dists = [Normal(2.0, 0.5), Normal(2.5, 0.8), Normal(1.8, 0.6)]
# Create weighted distributions with missing constructor weights
weighted_dists = weight(dists)
# Weights provided at observation time via joint observations
logpdf(weighted_dists, (values = y_obs, weights = [25, 10, 15]))
CensoredDistributions.weight
— Methodweight(
dist::Distributions.UnivariateDistribution,
weights::AbstractVector{<:Real}
) -> Any
Create a product distribution of weighted distributions, each with a different weight.
A Product
distribution of Weighted
distributions suitable for vectorized observations.
Arguments
dist
: The univariate distribution to be replicated and weighted for each observationweights
: Vector of weights to apply to each copy of the distribution
Examples
y_obs = [3.5, 4.2, 3.8] # Observed values
n_counts = [25, 10, 15] # Counts for each observation
d = Normal(2.0, 1.0)
weighted_dists = weight(d, n_counts)
# Weighted log-probability calculation
weighted_logpdf = logpdf(weighted_dists, y_obs)
# equivalent to: sum(n_counts .* logpdf.(d, y_obs))
See also
Weighted
: The underlying weighted distribution type
CensoredDistributions.weight
— Methodweight(
dist::Distributions.UnivariateDistribution,
w::Real
) -> CensoredDistributions.Weighted{D, T} where {D<:(Distributions.UnivariateDistribution), T<:Real}
Create a weighted distribution where the log-probability is scaled by w
.
A Weighted
distribution will contribute w * logpdf(dist, x)
to the log-probability when evaluating logpdf(weighted_dist, x)
.
Examples
# For aggregated count data
y_obs = 3.5 # Observed value
n_count = 25 # Number of times this value was observed
d = Normal(2.0, 1.0)
weighted_d = weight(d, n_count)
# Weighted log-probability calculation
weighted_logpdf = logpdf(weighted_d, y_obs)
# equivalent to: n_count * logpdf(d, y_obs)
CensoredDistributions.weight
— Methodweight(
dist::Distributions.UnivariateDistribution
) -> CensoredDistributions.Weighted{D, Missing} where D<:(Distributions.UnivariateDistribution)
Create a weighted distribution with missing constructor weight.
Useful for creating distributions where weights will be provided at observation time. Uses missing
as constructor weight, enabling observation weight to be used directly.
Examples
using CensoredDistributions, Distributions
d = Normal(2.0, 0.5)
# Create weighted distribution with missing constructor weight
weighted_dist = weight(d)
# Weight provided at observation time via joint observations
logpdf(weighted_dist, (value = 3.5, weight = 25))