epinowcast
Loading...
Searching...
No Matches
primary_censored_dist_analytical_cdf.stan
Go to the documentation of this file.
1
10int check_for_analytical(int dist_id, int primary_dist_id) {
11 if (dist_id == 2 && primary_dist_id == 1) return 1; // Gamma delay with Uniform primary
12 if (dist_id == 1 && primary_dist_id == 1) return 1; // Lognormal delay with Uniform primary
13 return 0; // No analytical solution for other combinations
14}
15
27real primary_censored_gamma_uniform_lcdf(data real d, real q, array[] real params, data real pwindow) {
28 real shape = params[1];
29 real rate = params[2];
30 real shape_1 = shape + 1;
31 real log_window = log(pwindow);
32
33 real log_F_T = gamma_lcdf(d | shape, rate);
34 real log_F_T_kp1 = gamma_lcdf(d | shape_1, rate);
35
36 real log_delta_F_T_kp1;
37 real log_delta_F_T_k;
38 real log_F_Splus;
39
40 if (q != 0) {
41 real log_F_T_q = gamma_lcdf(q | shape, rate);
42 real log_F_T_q_kp1 = gamma_lcdf(q | shape_1, rate);
43
44 // Ensure that the first argument is greater than the second
45 log_delta_F_T_kp1 = log_diff_exp(log_F_T_kp1, log_F_T_q_kp1);
46 log_delta_F_T_k = log_diff_exp(log_F_T, log_F_T_q);
47
48 log_F_Splus = log_diff_exp(
49 log_F_T,
50 log_diff_exp(
51 log(shape * inv(rate)) + log_delta_F_T_kp1,
52 log(d - pwindow) + log_delta_F_T_k
53 ) - log_window
54 );
55 } else {
56 log_delta_F_T_kp1 = log_F_T_kp1;
57 log_delta_F_T_k = log_F_T;
58
59 log_F_Splus = log_diff_exp(
60 log_F_T,
61 log_sum_exp(
62 log(shape * inv(rate)) + log_delta_F_T_kp1,
63 log(pwindow - d) + log_delta_F_T_k
64 ) - log_window
65 );
66 }
67
68 return log_F_Splus;
69}
70
82real primary_censored_lognormal_uniform_lcdf(data real d, real q, array[] real params, data real pwindow) {
83 real mu = params[1];
84 real sigma = params[2];
85 real mu_sigma2 = mu + square(sigma);
86 real log_window = log(pwindow);
87
88 real log_F_T = lognormal_lcdf(d | mu, sigma);
89 real log_F_T_mu_sigma2 = lognormal_lcdf(d | mu_sigma2, sigma);
90
91 real log_delta_F_T_mu_sigma;
92 real log_delta_F_T;
93 real log_F_Splus;
94
95 if (q != 0) {
96 real log_F_T_q = lognormal_lcdf(q | mu, sigma);
97 real log_F_T_q_mu_sigma2 = lognormal_lcdf(q | mu_sigma2, sigma);
98
99 // Ensure that the first argument is greater than the second
100 log_delta_F_T_mu_sigma = log_diff_exp(
101 log_F_T_mu_sigma2, log_F_T_q_mu_sigma2
102 );
103 log_delta_F_T = log_diff_exp(log_F_T, log_F_T_q);
104
105 log_F_Splus = log_diff_exp(
106 log_F_T,
107 log_diff_exp(
108 (mu + 0.5 * square(sigma)) + log_delta_F_T_mu_sigma,
109 log(d - pwindow) + log_delta_F_T
110 ) - log_window
111 );
112 } else {
113 log_delta_F_T_mu_sigma = log_F_T_mu_sigma2;
114 log_delta_F_T = log_F_T;
115
116 log_F_Splus = log_diff_exp(
117 log_F_T,
118 log_sum_exp(
119 (mu + 0.5 * square(sigma)) + log_delta_F_T_mu_sigma,
120 log(pwindow - d) + log_delta_F_T
121 ) - log_window
122 );
123 }
124
125 return log_F_Splus;
126}
127
141real primary_censored_dist_analytical_lcdf(data real d, int dist_id,
142 array[] real params,
143 data real pwindow, data real D,
144 int primary_dist_id,
145 array[] real primary_params) {
146 real result;
147 real log_cdf_D;
148
149 if (d <= 0) return negative_infinity();
150 if (d >= D) return 0;
151
152 real q = max({d - pwindow, 0});
153
154 if (dist_id == 2 && primary_dist_id == 1) {
155 // Gamma delay with Uniform primary
156 result = primary_censored_gamma_uniform_lcdf(d | q, params, pwindow);
157 } else if (dist_id == 1 && primary_dist_id == 1) {
158 // Lognormal delay with Uniform primary
159 result = primary_censored_lognormal_uniform_lcdf(d | q, params, pwindow);
160 } else {
161 // No analytical solution available
162 return negative_infinity();
163 }
164
165 if (!is_inf(D)) {
166 log_cdf_D = primary_censored_dist_lcdf(
167 D | dist_id, params, pwindow, positive_infinity(),
168 primary_dist_id, primary_params
169 );
170 result = result - log_cdf_D;
171 }
172
173 return result;
174}
175
189real primary_censored_dist_analytical_cdf(data real d, int dist_id,
190 array[] real params,
191 data real pwindow, data real D,
192 int primary_dist_id,
193 array[] real primary_params) {
194 return exp(primary_censored_dist_analytical_lcdf(d | dist_id, params, pwindow, D, primary_dist_id, primary_params));
195}
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_lognormal_uniform_lcdf(data real d, real q, array[] real params, data real pwindow)
int check_for_analytical(int dist_id, int primary_dist_id)
real primary_censored_dist_analytical_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_analytical_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_gamma_uniform_lcdf(data real d, real q, array[] real params, data real pwindow)