28 if (x < xmin || x > xmax) {
29 return negative_infinity();
32 return -log(xmax - xmin);
34 return log(abs(r)) + r * x -
35 log(abs(exp(r * xmax) - exp(r * xmin)));
38 if (dist_id == 2 && primary_id == 1)
return 1;
39 if (dist_id == 1 && primary_id == 1)
return 1;
40 if (dist_id == 3 && primary_id == 1)
return 1;
44 real shape = params[1];
45 real rate = params[2];
46 real log_window = log(pwindow);
48 real log_E = log(shape) - log(rate);
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);
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;
70 log_q_F_T_q = negative_infinity();
71 log_E_tF_T_q = negative_infinity();
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);
81 return log_diff_exp(log_A, log_B) - log_window;
85 real sigma = params[2];
86 real mu_sigma2 = mu + square(sigma);
87 real log_window = log(pwindow);
89 real log_E = mu + 0.5 * square(sigma);
91 real log_F_T_d = lognormal_lcdf(d | mu, sigma);
92 real log_tF_T_d = lognormal_lcdf(d | mu_sigma2, sigma);
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;
103 log_q_F_T_q = negative_infinity();
104 log_E_tF_T_q = negative_infinity();
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);
114 return log_diff_exp(log_A, log_B) - log_window;
117 real x = pow(t * inv(scale), shape);
118 real a = 1 + inv(shape);
119 return log(gamma_p(a, x)) + lgamma(a);
122 real shape = params[1];
123 real scale = params[2];
124 real log_window = log(pwindow);
125 real log_scale = log(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);
136 log_q_F_T_q = log(q) + weibull_lcdf(q | shape, scale);
139 log_q_F_T_q = negative_infinity();
140 log_E_tF_T_q = negative_infinity();
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);
150 return log_diff_exp(log_A, log_B) - log_window;
156 real q = max({d - pwindow, 0});
158 if (dist_id == 2 && primary_id == 1) {
160 }
else if (dist_id == 1 && primary_id == 1) {
162 }
else if (dist_id == 3 && primary_id == 1) {
165 return negative_infinity();
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;
176 d, dist_id, params, pwindow, primary_id
180 if (!is_inf(D) || L > 0) {
182 L, D, dist_id, params, pwindow, primary_id, primary_params
184 real log_cdf_L = bounds[1];
185 real log_cdf_D = bounds[2];
195 data real pwindow, data real L,
196 data real D,
int primary_id,
197 array[] real primary_params) {
201 if (dist_id == 1)
return 1;
202 if (dist_id == 2)
return 1;
203 if (dist_id == 3)
return 1;
204 if (dist_id == 4)
return 1;
205 if (dist_id == 9)
return 1;
206 if (dist_id == 13)
return 1;
207 if (dist_id == 16)
return 1;
208 if (dist_id == 19)
return 1;
209 if (dist_id == 21)
return 1;
210 if (dist_id == 22)
return 1;
213real
dist_lcdf(real delay, array[] real params,
int dist_id) {
215 return negative_infinity();
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);
239real
primary_lpdf(real x,
int primary_id, array[] real params, real xmin, real xmax) {
241 if (primary_id == 1)
return uniform_lpdf(x | xmin, xmax);
242 if (primary_id == 2)
return expgrowth_lpdf(x | xmin, xmax, params[1]);
244 reject(
"Invalid primary distribution identifier");
247 array[] real x_r, array[]
int x_i) {
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];
256 array[dist_params_len] real params;
257 if (dist_params_len) {
258 params = theta[1:dist_params_len];
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];
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);
269 return rep_vector(exp(log_cdf + log_primary_pdf), 1);
273 return log_diff_exp(log_cdf_D, log_cdf_L);
279 real log_normalizer, real L) {
281 return log_diff_exp(log_cdf, log_cdf_L) - log_normalizer;
283 return log_cdf - log_normalizer;
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
301 result[1] = negative_infinity();
304 L | dist_id, params, pwindow,
306 positive_infinity(), primary_id, primary_params
315 D | dist_id, params, pwindow,
317 positive_infinity(), primary_id, primary_params
324 data real pwindow, data real L, data real D,
326 array[] real primary_params) {
340 d | dist_id, params, pwindow, L, D, primary_id, primary_params
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};
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];
359 if (!is_inf(D) || !is_inf(L)) {
360 real log_result = log(result);
362 L, D, dist_id, params, pwindow, primary_id, primary_params
364 real log_cdf_L = bounds[1];
365 real log_cdf_D = bounds[2];
369 log_result, log_cdf_L, log_normalizer, L
371 result = exp(log_result);
378 data real pwindow, data real L, data real D,
380 array[] real primary_params) {
384 return negative_infinity();
396 d | dist_id, params, pwindow,
398 positive_infinity(), primary_id, primary_params
403 d | dist_id, params, pwindow,
405 positive_infinity(), primary_id, primary_params
410 if (!is_inf(D) || !is_inf(L)) {
412 L, D, dist_id, params, pwindow, primary_id, primary_params
414 real log_cdf_L = bounds[1];
415 real log_cdf_D = bounds[2];
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
429 int upper_interval = max_delay + 1;
430 vector[upper_interval] log_pmfs;
431 vector[upper_interval] log_cdfs;
435 if (D < upper_interval) {
436 reject(
"D must be at least max_delay + 1");
445 int start_idx = (!is_inf(L) && L > 0) ? max(1, to_int(floor(L))) : 1;
446 for (d in start_idx:upper_interval) {
448 d | dist_id, params, pwindow,
450 positive_infinity(), primary_id, primary_params
458 log_cdf_L = negative_infinity();
459 }
else if (L >= 1 && L <= upper_interval && floor(L) == L) {
461 log_cdf_L = log_cdfs[to_int(L)];
465 L | dist_id, params, pwindow,
467 positive_infinity(), primary_id, primary_params
473 if (D > upper_interval) {
478 D | dist_id, params, pwindow,
480 positive_infinity(), primary_id, primary_params
484 log_cdf_D = log_cdfs[upper_interval];
490 for (d in 1:upper_interval) {
493 log_pmfs[d] = negative_infinity();
494 }
else if (d - 1 < L) {
496 log_pmfs[d] = log_diff_exp(log_cdfs[d], log_cdf_L) - log_normalizer;
500 log_pmfs[d] = log_cdfs[d] - log_normalizer;
505 0.0 | dist_id, params, pwindow,
506 negative_infinity(), positive_infinity(),
507 primary_id, primary_params
509 log_pmfs[d] = log_diff_exp(log_cdfs[d], log_cdf_0) - log_normalizer;
512 log_pmfs[d] = log_diff_exp(log_cdfs[d], log_cdfs[d-1]) - log_normalizer;
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)