primarycensored
Loading...
Searching...
No Matches
primarycensored.stan
Go to the documentation of this file.
1
19real primarycensored_cdf(data real d, int dist_id, array[] real params,
20 data real pwindow, data real D,
21 int primary_id,
22 array[] real primary_params) {
23 real result;
24 if (d <= 0) {
25 return 0;
26 }
27
28 if (d >= D) {
29 return 1;
30 }
31
32 // Check if an analytical solution exists
33 if (check_for_analytical(dist_id, primary_id)) {
34 // Use analytical solution
36 d | dist_id, params, pwindow, D, primary_id, primary_params
37 );
38 } else {
39 // Use numerical integration for other cases
40 real lower_bound = max({d - pwindow, 1e-6});
41 int n_params = num_elements(params);
42 int n_primary_params = num_elements(primary_params);
43 array[n_params + n_primary_params] real theta = append_array(params, primary_params);
44 array[4] int ids = {dist_id, primary_id, n_params, n_primary_params};
45
46 vector[1] y0 = rep_vector(0.0, 1);
47 result = ode_rk45(primarycensored_ode, y0, lower_bound, {d}, theta, {d, pwindow}, ids)[1, 1];
48
49 if (!is_inf(D)) {
50 real log_cdf_D = primarycensored_lcdf(
51 D | dist_id, params, pwindow, positive_infinity(), primary_id,primary_params
52 );
53 result = exp(log(result) - log_cdf_D);
54 }
55 }
56
57 return result;
58}
59
88real primarycensored_lcdf(data real d, int dist_id, array[] real params,
89 data real pwindow, data real D,
90 int primary_id,
91 array[] real primary_params) {
92 real result;
93
94 if (d <= 0) {
95 return negative_infinity();
96 }
97
98 if (d >= D) {
99 return 0;
100 }
101
102 // Check if an analytical solution exists
103 if (check_for_analytical(dist_id, primary_id)) {
105 d | dist_id, params, pwindow, positive_infinity(), primary_id, primary_params
106 );
107 } else {
108 // Use numerical integration
109 result = log(primarycensored_cdf(
110 d | dist_id, params, pwindow, positive_infinity(), primary_id, primary_params
111 ));
112 }
113
114 // Handle truncation
115 if (!is_inf(D)) {
116 real log_cdf_D = primarycensored_lcdf(
117 D | dist_id, params, pwindow, positive_infinity(), primary_id, primary_params
118 );
119 result = result - log_cdf_D;
120 }
121
122 return result;
123}
124
155real primarycensored_lpmf(data int d, int dist_id, array[] real params,
156 data real pwindow, data real d_upper,
157 data real D, int primary_id,
158 array[] real primary_params) {
159 if (d_upper > D) {
160 reject("Upper truncation point is greater than D. It is ", d_upper,
161 " and D is ", D, ". Resolve this by increasing D to be greater or equal to d + swindow or decreasing swindow.");
162 }
163 if (d_upper <= d) {
164 reject("Upper truncation point is less than or equal to d. It is ", d_upper,
165 " and d is ", d, ". Resolve this by increasing d to be less than d_upper.");
166 }
167 real log_cdf_upper = primarycensored_lcdf(
168 d_upper | dist_id, params, pwindow, positive_infinity(), primary_id, primary_params
169 );
170 real log_cdf_lower = primarycensored_lcdf(
171 d | dist_id, params, pwindow, positive_infinity(), primary_id, primary_params
172 );
173 if (!is_inf(D)) {
174 real log_cdf_D;
175
176 if (d_upper == D) {
177 log_cdf_D = log_cdf_upper;
178 } else {
179 log_cdf_D = primarycensored_lcdf(
180 D | dist_id, params, pwindow, positive_infinity(), primary_id, primary_params
181 );
182 }
183 return log_diff_exp(log_cdf_upper, log_cdf_lower) - log_cdf_D;
184 } else {
185 return log_diff_exp(log_cdf_upper, log_cdf_lower);
186 }
187}
188
218real primarycensored_pmf(data int d, int dist_id, array[] real params,
219 data real pwindow, data real d_upper,
220 data real D, int primary_id,
221 array[] real primary_params) {
222 return exp(
224 d | dist_id, params, pwindow, d_upper, D, primary_id, primary_params
225 )
226 );
227}
228
269 int max_delay, data real D, int dist_id,
270 array[] real params, data real pwindow,
271 int primary_id, array[] real primary_params
272) {
273
274 int upper_interval = max_delay + 1;
275 vector[upper_interval] log_pmfs;
276 vector[upper_interval] log_cdfs;
277 real log_normalizer;
278
279 // Check if D is at least max_delay + 1
280 if (D < upper_interval) {
281 reject("D must be at least max_delay + 1");
282 }
283
284 // Compute log CDFs
285 for (d in 1:upper_interval) {
286 log_cdfs[d] = primarycensored_lcdf(
287 d | dist_id, params, pwindow, positive_infinity(), primary_id,
288 primary_params
289 );
290 }
291
292 // Compute log normalizer using upper_interval
293 if (D > upper_interval) {
294 if (is_inf(D)) {
295 log_normalizer = 0; // No normalization needed for infinite D
296 } else {
297 log_normalizer = primarycensored_lcdf(
298 D | dist_id, params, pwindow, positive_infinity(),
299 primary_id, primary_params
300 );
301 }
302 } else {
303 log_normalizer = log_cdfs[upper_interval];
304 }
305
306 // Compute log PMFs
307 log_pmfs[1] = log_cdfs[1] - log_normalizer;
308 for (d in 2:upper_interval) {
309 log_pmfs[d] = log_diff_exp(log_cdfs[d], log_cdfs[d-1]) - log_normalizer;
310 }
311
312 return log_pmfs;
313}
314
352 int max_delay, data real D, int dist_id,
353 array[] real params, data real pwindow,
354 int primary_id,
355 array[] real primary_params
356) {
357 return exp(
359 max_delay, D, dist_id, params, pwindow, primary_id, primary_params
360 )
361 );
362}
int check_for_analytical(int dist_id, int primary_id)
vector primarycensored_ode(real t, vector y, array[] real theta, array[] real x_r, array[] int x_i)
real primarycensored_cdf(data real d, int dist_id, array[] real params, data real pwindow, data real D, int primary_id, array[] real primary_params)
real primarycensored_lpmf(data int d, int dist_id, array[] real params, data real pwindow, data real d_upper, data real D, int primary_id, array[] real primary_params)
real primarycensored_lcdf(data real d, int dist_id, array[] real params, data real pwindow, data real D, int primary_id, array[] real primary_params)
real primarycensored_pmf(data int d, int dist_id, array[] real params, data real pwindow, data real d_upper, data real D, int primary_id, array[] real primary_params)
vector primarycensored_sone_pmf_vectorized(int max_delay, data real D, int dist_id, array[] real params, data real pwindow, int primary_id, array[] real primary_params)
vector primarycensored_sone_lpmf_vectorized(int max_delay, data real D, int dist_id, array[] real params, data real pwindow, 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 D, int primary_id, array[] real primary_params)
real primarycensored_analytical_lcdf(data real d, int dist_id, array[] real params, data real pwindow, data real D, int primary_id, array[] real primary_params)