epinowcast
Loading...
Searching...
No Matches
primary_censored_dist.stan File Reference

Go to the source code of this file.

Functions

real dist_lcdf (real delay, array[] real params, int dist_id)
 
real primary_dist_lpdf (real x, int primary_dist_id, array[] real params, real min, real max)
 
vector primary_censored_ode (real t, vector y, array[] real theta, array[] real x_r, array[] int x_i)
 
real primary_censored_dist_cdf (data real d, int dist_id, array[] real params, data real pwindow, data real D, int primary_dist_id, array[] real primary_params)
 
real primary_censored_dist_lcdf (data real d, int dist_id, array[] real params, data real pwindow, data real D, int primary_dist_id, array[] real primary_params)
 
real primary_censored_dist_lpmf (data int d, int dist_id, array[] real params, data real pwindow, data real d_upper, data real D, int primary_dist_id, array[] real primary_params)
 
real primary_censored_dist_pmf (data int d, int dist_id, array[] real params, data real pwindow, data real d_upper, data real D, int primary_dist_id, array[] real primary_params)
 
vector primary_censored_sone_lpmf_vectorized (int max_delay, data real D, int dist_id, array[] real params, data real pwindow, int primary_dist_id, array[] real primary_params)
 
vector primary_censored_sone_pmf_vectorized (int max_delay, data real D, int dist_id, array[] real params, data real pwindow, int primary_dist_id, array[] real primary_params)
 

Function Documentation

◆ dist_lcdf()

real dist_lcdf ( real  delay,
array[]real  params,
int  dist_id 
)

Primary event censored distribution functions Compute the log CDF of the delay distribution

Parameters
delayTime delay
paramsDistribution parameters
dist_idDistribution identifier 1: Lognormal, 2: Gamma, 3: Normal, 4: Exponential, 5: Weibull, 6: Beta, 7: Cauchy, 8: Chi-square, 9: Inverse Chi-square, 10: Double Exponential, 11: Inverse Gamma, 12: Logistic, 13: Pareto, 14: Scaled Inverse Chi-square, 15: Student's t, 16: Uniform, 17: von Mises
Returns
Log CDF of the delay distribution
// Example: Lognormal distribution
real delay = 5.0;
array[2] real params = {0.0, 1.0}; // mean and standard deviation on log scale
int dist_id = 1; // Lognormal
real log_cdf = dist_lcdf(delay, params, dist_id);
real dist_lcdf(real delay, array[] real params, int dist_id)

Definition at line 27 of file primary_censored_dist.stan.

◆ primary_censored_dist_cdf()

real primary_censored_dist_cdf ( data real  d,
int  dist_id,
array[]real  params,
data real  pwindow,
data real  D,
int  primary_dist_id,
array[]real  primary_params 
)

Compute the primary event censored CDF for a single delay

Parameters
dDelay
dist_idDistribution identifier
paramsArray of distribution parameters
pwindowPrimary event window
DMaximum delay (truncation point)
primary_dist_idPrimary distribution identifier
primary_paramsPrimary distribution parameters
Returns
Primary event censored CDF, normalized by D if finite (truncation adjustment)

Definition at line 130 of file primary_censored_dist.stan.

◆ primary_censored_dist_lcdf()

real primary_censored_dist_lcdf ( data real  d,
int  dist_id,
array[]real  params,
data real  pwindow,
data real  D,
int  primary_dist_id,
array[]real  primary_params 
)

Compute the primary event censored log CDF for a single delay

Parameters
dDelay
dist_idDistribution identifier
paramsArray of distribution parameters
pwindowPrimary event window
DMaximum delay (truncation point)
primary_dist_idPrimary distribution identifier
primary_paramsPrimary distribution parameters
Returns
Primary event censored log CDF, normalized by D if finite (truncation adjustment)
// Example: Weibull delay distribution with uniform primary distribution
real d = 3.0;
int dist_id = 5; // Weibull
array[2] real params = {2.0, 1.5}; // shape and scale
real pwindow = 1.0;
real D = positive_infinity();
int primary_dist_id = 1; // Uniform
array[0] real primary_params = {};
d, dist_id, params, pwindow, D, primary_dist_id, primary_params
);
real primary_censored_dist_lcdf(data real d, int dist_id, array[] real params, data real pwindow, data real D, int primary_dist_id, array[] real primary_params)

Definition at line 196 of file primary_censored_dist.stan.

◆ primary_censored_dist_lpmf()

real primary_censored_dist_lpmf ( data int  d,
int  dist_id,
array[]real  params,
data real  pwindow,
data real  d_upper,
data real  D,
int  primary_dist_id,
array[]real  primary_params 
)

Compute the primary event censored log PMF for a single delay

Parameters
dDelay (integer)
dist_idDistribution identifier
paramsArray of distribution parameters
pwindowPrimary event window
d_upperUpper bound for the delay interval
DMaximum delay (truncation point)
primary_dist_idPrimary distribution identifier
primary_paramsPrimary distribution parameters
Returns
Primary event censored log PMF, normalized by D if finite (truncation adjustment)
// Example: Weibull delay distribution with uniform primary distribution
int d = 3;
int dist_id = 5; // Weibull
array[2] real params = {2.0, 1.5}; // shape and scale
real pwindow = 1.0;
real d_upper = 4.0;
real D = positive_infinity();
int primary_dist_id = 1; // Uniform
array[0] real primary_params = {};
d, dist_id, params, pwindow, d_upper, D, primary_dist_id, primary_params
);
real primary_censored_dist_lpmf(data int d, int dist_id, array[] real params, data real pwindow, data real d_upper, data real D, int primary_dist_id, array[] real primary_params)

Definition at line 262 of file primary_censored_dist.stan.

◆ primary_censored_dist_pmf()

real primary_censored_dist_pmf ( data int  d,
int  dist_id,
array[]real  params,
data real  pwindow,
data real  d_upper,
data real  D,
int  primary_dist_id,
array[]real  primary_params 
)

Compute the primary event censored PMF for a single delay

Parameters
dDelay (integer)
dist_idDistribution identifier
paramsArray of distribution parameters
pwindowPrimary event window
d_upperUpper bound for the delay interval
DMaximum delay (truncation point)
primary_dist_idPrimary distribution identifier
primary_paramsPrimary distribution parameters
Returns
Primary event censored PMF, normalized by D if finite (truncation adjustment)
// Example: Weibull delay distribution with uniform primary distribution
int d = 3;
real d = 3.0;
int dist_id = 5; // Weibull
array[2] real params = {2.0, 1.5}; // shape and scale
real pwindow = 1.0;
real swindow = 0.1;
real D = positive_infinity();
int primary_dist_id = 1; // Uniform
array[0] real primary_params = {};
real pmf = primary_censored_dist_pmf(d, dist_id, params, pwindow, swindow, D, primary_dist_id, primary_params);
real primary_censored_dist_pmf(data int d, int dist_id, array[] real params, data real pwindow, data real d_upper, data real D, int primary_dist_id, array[] real primary_params)

Definition at line 324 of file primary_censored_dist.stan.

◆ primary_censored_ode()

vector primary_censored_ode ( real  t,
vector  y,
array[]real  theta,
array[]real  x_r,
array[]int  x_i 
)

ODE system for the primary censored distribution

Parameters
tTime
yState variables
thetaParameters
x_rReal data
x_iInteger data
Returns
Derivatives of the state variables

Definition at line 91 of file primary_censored_dist.stan.

◆ primary_censored_sone_lpmf_vectorized()

vector primary_censored_sone_lpmf_vectorized ( int  max_delay,
data real  D,
int  dist_id,
array[]real  params,
data real  pwindow,
int  primary_dist_id,
array[]real  primary_params 
)

Compute the primary event censored log PMF for integer delays up to max_delay

Parameters
max_delayMaximum delay to compute PMF for
DMaximum delay (truncation point), must be at least max_delay + 1
dist_idDistribution identifier
paramsArray of distribution parameters
pwindowPrimary event window
primary_dist_idPrimary distribution identifier
primary_paramsPrimary distribution parameters
Returns
Vector of primary event censored log PMFs for delays [0, 1] to [max_delay, max_delay + 1].

This function differs from primary_censored_dist_lpmf in that it:

  1. Computes PMFs for all integer delays from [0, 1] to [max_delay, max_delay + 1] in one call.
  2. Assumes integer delays (swindow = 1)
  3. Is more computationally efficient for multiple delay calculation as it reduces the number of integration calls.
// Example: Weibull delay distribution with uniform primary distribution
int max_delay = 10;
real D = 15.0;
int dist_id = 5; // Weibull
array[2] real params = {2.0, 1.5}; // shape and scale
real pwindow = 7.0;
int primary_dist_id = 1; // Uniform
array[0] real primary_params = {};
vector[max_delay] log_pmf =
max_delay, D, dist_id, params, pwindow, primary_dist_id,
primary_params
);
vector primary_censored_sone_lpmf_vectorized(int max_delay, data real D, int dist_id, array[] real params, data real pwindow, int primary_dist_id, array[] real primary_params)

Definition at line 373 of file primary_censored_dist.stan.

◆ primary_censored_sone_pmf_vectorized()

vector primary_censored_sone_pmf_vectorized ( int  max_delay,
data real  D,
int  dist_id,
array[]real  params,
data real  pwindow,
int  primary_dist_id,
array[]real  primary_params 
)

Compute the primary event censored PMF for integer delays up to max_delay

Parameters
max_delayMaximum delay to compute PMF for
DMaximum delay (truncation point), must be at least max_delay + 1
dist_idDistribution identifier
paramsArray of distribution parameters
pwindowPrimary event window
primary_dist_idPrimary distribution identifier
primary_paramsPrimary distribution parameters
Returns
Vector of primary event censored PMFs for integer delays 1 to max_delay

This function differs from primary_censored_dist_pmf in that it:

  1. Computes PMFs for all integer delays from [0, 1] to [max_delay, max_delay + 1] in one call.
  2. Assumes integer delays (swindow = 1)
  3. Is more computationally efficient for multiple delay calculations
// Example: Weibull delay distribution with uniform primary distribution
int max_delay = 10;
real D = 15.0;
int dist_id = 5; // Weibull
array[2] real params = {2.0, 1.5}; // shape and scale
real pwindow = 7.0;
int primary_dist_id = 1; // Uniform
array[0] real primary_params = {};
vector[max_delay] pmf =
max_delay, D, dist_id, params, pwindow, primary_dist_id, primary_params
);
vector primary_censored_sone_pmf_vectorized(int max_delay, data real D, int dist_id, array[] real params, data real pwindow, int primary_dist_id, array[] real primary_params)

Definition at line 455 of file primary_censored_dist.stan.

◆ primary_dist_lpdf()

real primary_dist_lpdf ( real  x,
int  primary_dist_id,
array[]real  params,
real  min,
real  max 
)

Compute the log PDF of the primary distribution

Parameters
xValue
primary_dist_idPrimary distribution identifier
paramsDistribution parameters
minMinimum value
maxMaximum value
Returns
Log PDF of the primary distribution
// Example: Uniform distribution
real x = 0.5;
int primary_dist_id = 1; // Uniform
array[0] real params = {}; // No additional parameters for uniform
real min = 0;
real max = 1;
real log_pdf = primary_dist_lpdf(x, primary_dist_id, params, min, max);
real primary_dist_lpdf(real x, int primary_dist_id, array[] real params, real min, real max)

Definition at line 72 of file primary_censored_dist.stan.