epinowcast
Loading...
Searching...
No Matches
primarycensored.stan
Go to the documentation of this file.
1
26// Stan functions from primarycensored version 1.5.0
27real expgrowth_lpdf(real x, real xmin, real xmax, real r) {
28 if (x < xmin || x > xmax) {
29 return negative_infinity();
30 }
31 if (abs(r) < 1e-10) {
32 return -log(xmax - xmin);
33 }
34 return log(abs(r)) + r * x -
35 log(abs(exp(r * xmax) - exp(r * xmin)));
36}
37int check_for_analytical(int dist_id, int primary_id) {
38 if (dist_id == 2 && primary_id == 1) return 1; // Gamma delay with Uniform primary
39 if (dist_id == 1 && primary_id == 1) return 1; // Lognormal delay with Uniform primary
40 if (dist_id == 3 && primary_id == 1) return 1; // Weibull delay with Uniform primary
41 return 0; // No analytical solution for other combinations
42}
43real primarycensored_gamma_uniform_lcdf(data real d, real q, array[] real params, data real pwindow) {
44 real shape = params[1];
45 real rate = params[2];
46 real log_window = log(pwindow);
47 // log E where E = k * theta = shape / rate is the mean of the delay
48 real log_E = log(shape) - log(rate);
49
50 // F_T(d; k) and the recursion to F_T(d; k+1):
51 // P(k+1, y) = P(k, y) - y^k e^{-y} / Gamma(k+1), with y = rate * d
52 real log_F_T_d_k = gamma_lcdf(d | shape, rate);
53 real gamma_kp1_pdf_log_d
54 = shape * log(rate * d) - rate * d - lgamma(shape + 1);
55 real log_F_T_d_kp1 = log_diff_exp(log_F_T_d_k, gamma_kp1_pdf_log_d);
56
57 // q-dependent terms. Final algebra is unified; only a guard to avoid
58 // log_diff_exp(-inf, -inf) and log(0) when q == 0 (q is data, so autodiff
59 // is unaffected by this branch).
60 real log_q_F_T_q; // log(q * F_T(q; k))
61 real log_E_tF_T_q; // log(E * F_T(q; k+1))
62 if (q > 0) {
63 real log_F_T_q_k = gamma_lcdf(q | shape, rate);
64 real gamma_kp1_pdf_log_q
65 = shape * log(rate * q) - rate * q - lgamma(shape + 1);
66 real log_F_T_q_kp1 = log_diff_exp(log_F_T_q_k, gamma_kp1_pdf_log_q);
67 log_q_F_T_q = log(q) + log_F_T_q_k;
68 log_E_tF_T_q = log_E + log_F_T_q_kp1;
69 } else {
70 log_q_F_T_q = negative_infinity();
71 log_E_tF_T_q = negative_infinity();
72 }
73
74 // Unified form: F_{S+}(d) = (A - B) / w_P with A, B sums of positives:
75 // A = d * F_T(d; k) + E * F_T(q; k+1)
76 // B = q * F_T(q; k) + E * F_T(d; k+1)
77 // Ordering A >= B is guaranteed by F_{S+}(d) >= 0.
78 real log_A = log_sum_exp(log(d) + log_F_T_d_k, log_E_tF_T_q);
79 real log_B = log_sum_exp(log_q_F_T_q, log_E + log_F_T_d_kp1);
80
81 return log_diff_exp(log_A, log_B) - log_window;
82}
83real primarycensored_lognormal_uniform_lcdf(data real d, real q, array[] real params, data real pwindow) {
84 real mu = params[1];
85 real sigma = params[2];
86 real mu_sigma2 = mu + square(sigma);
87 real log_window = log(pwindow);
88 // log E where E = exp(mu + sigma^2/2) is the mean of the delay
89 real log_E = mu + 0.5 * square(sigma);
90
91 real log_F_T_d = lognormal_lcdf(d | mu, sigma);
92 real log_tF_T_d = lognormal_lcdf(d | mu_sigma2, sigma);
93
94 // q-dependent terms (guard only to avoid log(0); final algebra is unified).
95 real log_q_F_T_q; // log(q * F_T(q))
96 real log_E_tF_T_q; // log(E * tilde F_T(q))
97 if (q > 0) {
98 real log_F_T_q = lognormal_lcdf(q | mu, sigma);
99 real log_tF_T_q = lognormal_lcdf(q | mu_sigma2, sigma);
100 log_q_F_T_q = log(q) + log_F_T_q;
101 log_E_tF_T_q = log_E + log_tF_T_q;
102 } else {
103 log_q_F_T_q = negative_infinity();
104 log_E_tF_T_q = negative_infinity();
105 }
106
107 // Unified form: F_{S+}(d) = (A - B) / w_P with
108 // A = d * F_T(d) + E * tilde F_T(q)
109 // B = q * F_T(q) + E * tilde F_T(d)
110 // Ordering A >= B is guaranteed by F_{S+}(d) >= 0.
111 real log_A = log_sum_exp(log(d) + log_F_T_d, log_E_tF_T_q);
112 real log_B = log_sum_exp(log_q_F_T_q, log_E + log_tF_T_d);
113
114 return log_diff_exp(log_A, log_B) - log_window;
115}
116real log_weibull_g(real t, real shape, real scale) {
117 real x = pow(t * inv(scale), shape);
118 real a = 1 + inv(shape);
119 return log(gamma_p(a, x)) + lgamma(a);
120}
121real primarycensored_weibull_uniform_lcdf(data real d, real q, array[] real params, data real pwindow) {
122 real shape = params[1];
123 real scale = params[2];
124 real log_window = log(pwindow);
125 real log_scale = log(scale);
126
127 // For Weibull: E = scale (lambda) and tilde F_T(t) = g(t; lambda, k), so
128 // log(E * tilde F_T(t)) = log(scale) + log_weibull_g(t, shape, scale).
129 real log_F_T_d = weibull_lcdf(d | shape, scale);
130 real log_E_tF_T_d = log_scale + log_weibull_g(d, shape, scale);
131
132 // q-dependent terms (guard only to avoid log(0); final algebra is unified).
133 real log_q_F_T_q; // log(q * F_T(q))
134 real log_E_tF_T_q; // log(E * tilde F_T(q)) = log(scale * g(q; lambda, k))
135 if (q > 0) {
136 log_q_F_T_q = log(q) + weibull_lcdf(q | shape, scale);
137 log_E_tF_T_q = log_scale + log_weibull_g(q, shape, scale);
138 } else {
139 log_q_F_T_q = negative_infinity();
140 log_E_tF_T_q = negative_infinity();
141 }
142
143 // Unified form: F_{S+}(d) = (A - B) / w_P with
144 // A = d * F_T(d) + scale * g(q; lambda, k)
145 // B = q * F_T(q) + scale * g(d; lambda, k)
146 // Ordering A >= B is guaranteed by F_{S+}(d) >= 0.
147 real log_A = log_sum_exp(log(d) + log_F_T_d, log_E_tF_T_q);
148 real log_B = log_sum_exp(log_q_F_T_q, log_E_tF_T_d);
149
150 return log_diff_exp(log_A, log_B) - log_window;
151}
152real primarycensored_analytical_lcdf_raw(data real d, int dist_id,
153 array[] real params,
154 data real pwindow,
155 int primary_id) {
156 real q = max({d - pwindow, 0});
157
158 if (dist_id == 2 && primary_id == 1) {
159 return primarycensored_gamma_uniform_lcdf(d | q, params, pwindow);
160 } else if (dist_id == 1 && primary_id == 1) {
161 return primarycensored_lognormal_uniform_lcdf(d | q, params, pwindow);
162 } else if (dist_id == 3 && primary_id == 1) {
163 return primarycensored_weibull_uniform_lcdf(d | q, params, pwindow);
164 }
165 return negative_infinity();
166}
167real primarycensored_analytical_lcdf(data real d, int dist_id,
168 array[] real params,
169 data real pwindow, data real L,
170 data real D, int primary_id,
171 array[] real primary_params) {
172 if (d <= L) return negative_infinity();
173 if (d >= D) return 0;
174
176 d, dist_id, params, pwindow, primary_id
177 );
178
179 // Apply truncation normalization
180 if (!is_inf(D) || L > 0) {
181 vector[2] bounds = primarycensored_truncation_bounds(
182 L, D, dist_id, params, pwindow, primary_id, primary_params
183 );
184 real log_cdf_L = bounds[1];
185 real log_cdf_D = bounds[2];
186
187 real log_normalizer = primarycensored_log_normalizer(log_cdf_D, log_cdf_L, L);
188 result = primarycensored_apply_truncation(result, log_cdf_L, log_normalizer, L);
189 }
190
191 return result;
192}
193real primarycensored_analytical_cdf(data real d, int dist_id,
194 array[] real params,
195 data real pwindow, data real L,
196 data real D, int primary_id,
197 array[] real primary_params) {
198 return exp(primarycensored_analytical_lcdf(d | dist_id, params, pwindow, L, D, primary_id, primary_params));
199}
200int dist_has_positive_support(data int dist_id) {
201 if (dist_id == 1) return 1; // Lognormal
202 if (dist_id == 2) return 1; // Gamma
203 if (dist_id == 3) return 1; // Weibull
204 if (dist_id == 4) return 1; // Exponential
205 if (dist_id == 9) return 1; // Beta (support on [0, 1])
206 if (dist_id == 13) return 1; // Chi-square
207 if (dist_id == 16) return 1; // Inverse Gamma
208 if (dist_id == 19) return 1; // Inverse Chi-square
209 if (dist_id == 21) return 1; // Pareto
210 if (dist_id == 22) return 1; // Scaled inverse Chi-square
211 return 0;
212}
213real dist_lcdf(real delay, array[] real params, int dist_id) {
214 if (dist_has_positive_support(dist_id) && delay <= 0) {
215 return negative_infinity();
216 }
217
218 // IDs match pcd_distributions$stan_id in R
219 if (dist_id == 1) return lognormal_lcdf(delay | params[1], params[2]);
220 else if (dist_id == 2) return gamma_lcdf(delay | params[1], params[2]);
221 else if (dist_id == 3) return weibull_lcdf(delay | params[1], params[2]);
222 else if (dist_id == 4) return exponential_lcdf(delay | params[1]);
223 else if (dist_id == 9) return beta_lcdf(delay | params[1], params[2]);
224 else if (dist_id == 12) return cauchy_lcdf(delay | params[1], params[2]);
225 else if (dist_id == 13) return chi_square_lcdf(delay | params[1]);
226 else if (dist_id == 15) return gumbel_lcdf(delay | params[1], params[2]);
227 else if (dist_id == 16) return inv_gamma_lcdf(delay | params[1], params[2]);
228 else if (dist_id == 17) return logistic_lcdf(delay | params[1], params[2]);
229 else if (dist_id == 18) return normal_lcdf(delay | params[1], params[2]);
230 else if (dist_id == 19) return inv_chi_square_lcdf(delay | params[1]);
231 else if (dist_id == 20) return double_exponential_lcdf(delay | params[1], params[2]);
232 else if (dist_id == 21) return pareto_lcdf(delay | params[1], params[2]);
233 else if (dist_id == 22) return scaled_inv_chi_square_lcdf(delay | params[1], params[2]);
234 else if (dist_id == 23) return student_t_lcdf(delay | params[1], params[2], params[3]);
235 else if (dist_id == 24) return uniform_lcdf(delay | params[1], params[2]);
236 else if (dist_id == 25) return von_mises_lcdf(delay | params[1], params[2]);
237 else reject("Invalid distribution identifier: ", dist_id);
238}
239real primary_lpdf(real x, int primary_id, array[] real params, real xmin, real xmax) {
240 // Implement switch for different primary distributions
241 if (primary_id == 1) return uniform_lpdf(x | xmin, xmax);
242 if (primary_id == 2) return expgrowth_lpdf(x | xmin, xmax, params[1]);
243 // Add more primary distributions as needed
244 reject("Invalid primary distribution identifier");
245}
246vector primarycensored_ode(real t, vector y, array[] real theta,
247 array[] real x_r, array[] int x_i) {
248 real d = x_r[1];
249 int dist_id = x_i[1];
250 int primary_id = x_i[2];
251 real pwindow = x_r[2];
252 int dist_params_len = x_i[3];
253 int primary_params_len = x_i[4];
254
255 // Extract distribution parameters
256 array[dist_params_len] real params;
257 if (dist_params_len) {
258 params = theta[1:dist_params_len];
259 }
260 array[primary_params_len] real primary_params;
261 if (primary_params_len) {
262 int primary_loc = num_elements(theta);
263 primary_params = theta[primary_loc - primary_params_len + 1:primary_loc];
264 }
265
266 real log_cdf = dist_lcdf(t | params, dist_id);
267 real log_primary_pdf = primary_lpdf(d - t | primary_id, primary_params, 0, pwindow);
268
269 return rep_vector(exp(log_cdf + log_primary_pdf), 1);
270}
271real primarycensored_log_normalizer(real log_cdf_D, real log_cdf_L, real L) {
272 if (!is_inf(L)) {
273 return log_diff_exp(log_cdf_D, log_cdf_L);
274 } else {
275 return log_cdf_D;
276 }
277}
278real primarycensored_apply_truncation(real log_cdf, real log_cdf_L,
279 real log_normalizer, real L) {
280 if (!is_inf(L)) {
281 return log_diff_exp(log_cdf, log_cdf_L) - log_normalizer;
282 } else {
283 return log_cdf - log_normalizer;
284 }
285}
287 data real L, data real D,
288 data int dist_id, array[] real params, data real pwindow,
289 data int primary_id, array[] real primary_params
290) {
291 vector[2] result;
292 // Internal lower bound for the un-truncated distribution: 0 lets the
293 // `d <= L` early-exit in primarycensored_lcdf return -inf for delays below
294 // the natural support of positive-support distributions; -inf disables that
295 // short-circuit so distributions with support on the reals are integrated.
296 // Expression is inlined (rather than bound to a local) so Stan's data-flow
297 // checker recognises it as data-only.
298
299 // Get CDF at lower truncation point L
300 if (is_inf(L)) {
301 result[1] = negative_infinity();
302 } else {
303 result[1] = primarycensored_lcdf(
304 L | dist_id, params, pwindow,
305 dist_has_positive_support(dist_id) ? 0.0 : negative_infinity(),
306 positive_infinity(), primary_id, primary_params
307 );
308 }
309
310 // Get CDF at upper truncation point D
311 if (is_inf(D)) {
312 result[2] = 0;
313 } else {
314 result[2] = primarycensored_lcdf(
315 D | dist_id, params, pwindow,
316 dist_has_positive_support(dist_id) ? 0.0 : negative_infinity(),
317 positive_infinity(), primary_id, primary_params
318 );
319 }
320
321 return result;
322}
323real primarycensored_cdf(data real d, data int dist_id, array[] real params,
324 data real pwindow, data real L, data real D,
325 data int primary_id,
326 array[] real primary_params) {
327 real result;
328 if (d <= L) {
329 return 0;
330 }
331
332 if (d >= D) {
333 return 1;
334 }
335
336 // Check if an analytical solution exists
337 if (check_for_analytical(dist_id, primary_id)) {
338 // Use analytical solution
340 d | dist_id, params, pwindow, L, D, primary_id, primary_params
341 );
342 } else {
343 // Use numerical integration for other cases. The integration variable
344 // ranges over the primary-event time, so the natural lower bound is
345 // d - pwindow. For positive-support delays the integrand `F_delay(t)` is
346 // 0 for t <= 0, so an unclipped lower bound just adds a flat zero region
347 // for negative t. Distributions with support on the reals also accept the
348 // unclipped lower bound directly.
349 real lower_bound = d - pwindow;
350 int n_params = num_elements(params);
351 int n_primary_params = num_elements(primary_params);
352 array[n_params + n_primary_params] real theta = append_array(params, primary_params);
353 array[4] int ids = {dist_id, primary_id, n_params, n_primary_params};
354
355 vector[1] y0 = rep_vector(0.0, 1);
356 result = ode_rk45(primarycensored_ode, y0, lower_bound, {d}, theta, {d, pwindow}, ids)[1, 1];
357
358 // Apply truncation normalization on log scale for numerical stability
359 if (!is_inf(D) || !is_inf(L)) {
360 real log_result = log(result);
361 vector[2] bounds = primarycensored_truncation_bounds(
362 L, D, dist_id, params, pwindow, primary_id, primary_params
363 );
364 real log_cdf_L = bounds[1];
365 real log_cdf_D = bounds[2];
366
367 real log_normalizer = primarycensored_log_normalizer(log_cdf_D, log_cdf_L, L);
369 log_result, log_cdf_L, log_normalizer, L
370 );
371 result = exp(log_result);
372 }
373 }
374
375 return result;
376}
377real primarycensored_lcdf(data real d, data int dist_id, array[] real params,
378 data real pwindow, data real L, data real D,
379 data int primary_id,
380 array[] real primary_params) {
381 real result;
382
383 if (d <= L) {
384 return negative_infinity();
385 }
386
387 if (d >= D) {
388 return 0;
389 }
390
391 // Check if an analytical solution exists. The internal lower bound is 0 for
392 // positive-support delays (lets the d <= L early-exit return -inf for d <= 0)
393 // and -inf for distributions with support on the reals.
394 if (check_for_analytical(dist_id, primary_id)) {
396 d | dist_id, params, pwindow,
397 dist_has_positive_support(dist_id) ? 0.0 : negative_infinity(),
398 positive_infinity(), primary_id, primary_params
399 );
400 } else {
401 // Use numerical integration
402 result = log(primarycensored_cdf(
403 d | dist_id, params, pwindow,
404 dist_has_positive_support(dist_id) ? 0.0 : negative_infinity(),
405 positive_infinity(), primary_id, primary_params
406 ));
407 }
408
409 // Handle truncation normalization
410 if (!is_inf(D) || !is_inf(L)) {
411 vector[2] bounds = primarycensored_truncation_bounds(
412 L, D, dist_id, params, pwindow, primary_id, primary_params
413 );
414 real log_cdf_L = bounds[1];
415 real log_cdf_D = bounds[2];
416
417 real log_normalizer = primarycensored_log_normalizer(log_cdf_D, log_cdf_L, L);
418 result = primarycensored_apply_truncation(result, log_cdf_L, log_normalizer, L);
419 }
420
421 return result;
422}
424 data int max_delay, data real L, data real D, data int dist_id,
425 array[] real params, data real pwindow,
426 data int primary_id, array[] real primary_params
427) {
428
429 int upper_interval = max_delay + 1;
430 vector[upper_interval] log_pmfs;
431 vector[upper_interval] log_cdfs;
432 real log_normalizer;
433
434 // Check if D is at least max_delay + 1
435 if (D < upper_interval) {
436 reject("D must be at least max_delay + 1");
437 }
438
439 // Compute log CDFs (without truncation normalization). The internal lower
440 // bound below is 0 for positive-support delays and -inf otherwise; it is
441 // inlined rather than bound to a local so Stan's type checker treats it as
442 // data-only.
443 // Start from max(1, floor(L)) to avoid computing unused CDFs when L > 0;
444 // for L <= 0 (including -inf) start at 1 since F(d) = 0 for d <= 0.
445 int start_idx = (!is_inf(L) && L > 0) ? max(1, to_int(floor(L))) : 1;
446 for (d in start_idx:upper_interval) {
447 log_cdfs[d] = primarycensored_lcdf(
448 d | dist_id, params, pwindow,
449 dist_has_positive_support(dist_id) ? 0.0 : negative_infinity(),
450 positive_infinity(), primary_id, primary_params
451 );
452 }
453
454 // Get CDF at lower truncation point L
455 real log_cdf_L;
456 if (is_inf(L)) {
457 // No left truncation (L = -inf sentinel)
458 log_cdf_L = negative_infinity();
459 } else if (L >= 1 && L <= upper_interval && floor(L) == L) {
460 // L is a positive integer within computed range, reuse cached value
461 log_cdf_L = log_cdfs[to_int(L)];
462 } else {
463 // L is outside computed range or non-integer, compute directly
464 log_cdf_L = primarycensored_lcdf(
465 L | dist_id, params, pwindow,
466 dist_has_positive_support(dist_id) ? 0.0 : negative_infinity(),
467 positive_infinity(), primary_id, primary_params
468 );
469 }
470
471 // Compute log normalizer: log(F(D) - F(L))
472 real log_cdf_D;
473 if (D > upper_interval) {
474 if (is_inf(D)) {
475 log_cdf_D = 0; // log(1) = 0 for infinite D
476 } else {
477 log_cdf_D = primarycensored_lcdf(
478 D | dist_id, params, pwindow,
479 dist_has_positive_support(dist_id) ? 0.0 : negative_infinity(),
480 positive_infinity(), primary_id, primary_params
481 );
482 }
483 } else {
484 log_cdf_D = log_cdfs[upper_interval];
485 }
486
487 log_normalizer = primarycensored_log_normalizer(log_cdf_D, log_cdf_L, L);
488
489 // Compute log PMFs: log((F(d) - F(d-1)) / (F(D) - F(L)))
490 for (d in 1:upper_interval) {
491 if (d <= L) {
492 // Delay interval [d-1, d) is entirely at or below L
493 log_pmfs[d] = negative_infinity();
494 } else if (d - 1 < L) {
495 // L falls within interval [d-1, d), so compute mass in [L, d)
496 log_pmfs[d] = log_diff_exp(log_cdfs[d], log_cdf_L) - log_normalizer;
497 } else if (d == 1 && dist_has_positive_support(dist_id)) {
498 // First interval [0, 1) with L <= 0 and positive-support delay:
499 // F(0) = 0, so PMF = F(1) / normalizer
500 log_pmfs[d] = log_cdfs[d] - log_normalizer;
501 } else if (d == 1) {
502 // First interval [0, 1) with L <= 0 and real-support delay: F(0) is
503 // non-zero in general, so compute it explicitly.
504 real log_cdf_0 = primarycensored_lcdf(
505 0.0 | dist_id, params, pwindow,
506 negative_infinity(), positive_infinity(),
507 primary_id, primary_params
508 );
509 log_pmfs[d] = log_diff_exp(log_cdfs[d], log_cdf_0) - log_normalizer;
510 } else {
511 // Standard case: PMF = (F(d) - F(d-1)) / normalizer
512 log_pmfs[d] = log_diff_exp(log_cdfs[d], log_cdfs[d-1]) - log_normalizer;
513 }
514 }
515
516 return log_pmfs;
517}
real dist_lcdf(real delay, array[] real params, int dist_id)
vector primarycensored_sone_lpmf_vectorized(data int max_delay, data real L, data real D, data int dist_id, array[] real params, data real pwindow, data int primary_id, array[] real primary_params)
real primarycensored_apply_truncation(real log_cdf, real log_cdf_L, real log_normalizer, real L)
real primarycensored_cdf(data real d, data int dist_id, array[] real params, data real pwindow, data real L, data real D, data int primary_id, array[] real primary_params)
real primarycensored_log_normalizer(real log_cdf_D, real log_cdf_L, real L)
real log_weibull_g(real t, real shape, real scale)
vector primarycensored_ode(real t, vector y, array[] real theta, array[] real x_r, array[] int x_i)
real primarycensored_analytical_lcdf_raw(data real d, int dist_id, array[] real params, data real pwindow, int primary_id)
real primarycensored_weibull_uniform_lcdf(data real d, real q, array[] real params, data real pwindow)
real expgrowth_lpdf(real x, real xmin, real xmax, real r)
real primarycensored_lcdf(data real d, data int dist_id, array[] real params, data real pwindow, data real L, data real D, data int primary_id, array[] real primary_params)
real primarycensored_lognormal_uniform_lcdf(data real d, real q, array[] real params, data real pwindow)
real primarycensored_analytical_lcdf(data real d, int dist_id, array[] real params, data real pwindow, data real L, data real D, int primary_id, array[] real primary_params)
int dist_has_positive_support(data int dist_id)
real primarycensored_gamma_uniform_lcdf(data real d, real q, array[] real params, data real pwindow)
real primary_lpdf(real x, int primary_id, array[] real params, real xmin, real xmax)
vector primarycensored_truncation_bounds(data real L, data real D, data int dist_id, array[] real params, data real pwindow, data int primary_id, array[] real primary_params)
real primarycensored_analytical_cdf(data real d, int dist_id, array[] real params, data real pwindow, data real L, data real D, int primary_id, array[] real primary_params)
int check_for_analytical(int dist_id, int primary_id)