| 1 | /* Copyright 2005 Robert Kern (robert.kern@gmail.com) |
|---|
| 2 | * |
|---|
| 3 | * Permission is hereby granted, free of charge, to any person obtaining a |
|---|
| 4 | * copy of this software and associated documentation files (the |
|---|
| 5 | * "Software"), to deal in the Software without restriction, including |
|---|
| 6 | * without limitation the rights to use, copy, modify, merge, publish, |
|---|
| 7 | * distribute, sublicense, and/or sell copies of the Software, and to |
|---|
| 8 | * permit persons to whom the Software is furnished to do so, subject to |
|---|
| 9 | * the following conditions: |
|---|
| 10 | * |
|---|
| 11 | * The above copyright notice and this permission notice shall be included |
|---|
| 12 | * in all copies or substantial portions of the Software. |
|---|
| 13 | * |
|---|
| 14 | * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS |
|---|
| 15 | * OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF |
|---|
| 16 | * MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. |
|---|
| 17 | * IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY |
|---|
| 18 | * CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, |
|---|
| 19 | * TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE |
|---|
| 20 | * SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE. |
|---|
| 21 | */ |
|---|
| 22 | |
|---|
| 23 | /* The implementations of rk_hypergeometric_hyp(), rk_hypergeometric_hrua(), |
|---|
| 24 | * and rk_triangular() were adapted from Ivan Frohne's rv.py which has this |
|---|
| 25 | * license: |
|---|
| 26 | * |
|---|
| 27 | * Copyright 1998 by Ivan Frohne; Wasilla, Alaska, U.S.A. |
|---|
| 28 | * All Rights Reserved |
|---|
| 29 | * |
|---|
| 30 | * Permission to use, copy, modify and distribute this software and its |
|---|
| 31 | * documentation for any purpose, free of charge, is granted subject to the |
|---|
| 32 | * following conditions: |
|---|
| 33 | * The above copyright notice and this permission notice shall be included in |
|---|
| 34 | * all copies or substantial portions of the software. |
|---|
| 35 | * |
|---|
| 36 | * THE SOFTWARE AND DOCUMENTATION IS PROVIDED WITHOUT WARRANTY OF ANY KIND, |
|---|
| 37 | * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO MERCHANTABILITY, FITNESS |
|---|
| 38 | * FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHOR |
|---|
| 39 | * OR COPYRIGHT HOLDER BE LIABLE FOR ANY CLAIM OR DAMAGES IN A CONTRACT |
|---|
| 40 | * ACTION, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE |
|---|
| 41 | * SOFTWARE OR ITS DOCUMENTATION. |
|---|
| 42 | */ |
|---|
| 43 | |
|---|
| 44 | #include <math.h> |
|---|
| 45 | #include "distributions.h" |
|---|
| 46 | #include <stdio.h> |
|---|
| 47 | |
|---|
| 48 | #ifndef min |
|---|
| 49 | #define min(x,y) ((x<y)?x:y) |
|---|
| 50 | #define max(x,y) ((x>y)?x:y) |
|---|
| 51 | #endif |
|---|
| 52 | |
|---|
| 53 | #ifndef M_PI |
|---|
| 54 | #define M_PI 3.14159265358979323846264338328 |
|---|
| 55 | #endif |
|---|
| 56 | /* log-gamma function to support some of these distributions. The |
|---|
| 57 | * algorithm comes from SPECFUN by Shanjie Zhang and Jianming Jin and their |
|---|
| 58 | * book "Computation of Special Functions", 1996, John Wiley & Sons, Inc. |
|---|
| 59 | */ |
|---|
| 60 | extern double loggam(double x); |
|---|
| 61 | double loggam(double x) |
|---|
| 62 | { |
|---|
| 63 | double x0, x2, xp, gl, gl0; |
|---|
| 64 | long k, n; |
|---|
| 65 | |
|---|
| 66 | static double a[10] = {8.333333333333333e-02,-2.777777777777778e-03, |
|---|
| 67 | 7.936507936507937e-04,-5.952380952380952e-04, |
|---|
| 68 | 8.417508417508418e-04,-1.917526917526918e-03, |
|---|
| 69 | 6.410256410256410e-03,-2.955065359477124e-02, |
|---|
| 70 | 1.796443723688307e-01,-1.39243221690590e+00}; |
|---|
| 71 | x0 = x; |
|---|
| 72 | n = 0; |
|---|
| 73 | if ((x == 1.0) || (x == 2.0)) |
|---|
| 74 | { |
|---|
| 75 | return 0.0; |
|---|
| 76 | } |
|---|
| 77 | else if (x <= 7.0) |
|---|
| 78 | { |
|---|
| 79 | n = (long)(7 - x); |
|---|
| 80 | x0 = x + n; |
|---|
| 81 | } |
|---|
| 82 | x2 = 1.0/(x0*x0); |
|---|
| 83 | xp = 2*M_PI; |
|---|
| 84 | gl0 = a[9]; |
|---|
| 85 | for (k=8; k>=0; k--) |
|---|
| 86 | { |
|---|
| 87 | gl0 *= x2; |
|---|
| 88 | gl0 += a[k]; |
|---|
| 89 | } |
|---|
| 90 | gl = gl0/x0 + 0.5*log(xp) + (x0-0.5)*log(x0) - x0; |
|---|
| 91 | if (x <= 7.0) |
|---|
| 92 | { |
|---|
| 93 | for (k=1; k<=n; k++) |
|---|
| 94 | { |
|---|
| 95 | gl -= log(x0-1.0); |
|---|
| 96 | x0 -= 1.0; |
|---|
| 97 | } |
|---|
| 98 | } |
|---|
| 99 | return gl; |
|---|
| 100 | } |
|---|
| 101 | |
|---|
| 102 | double rk_normal(rk_state *state, double loc, double scale) |
|---|
| 103 | { |
|---|
| 104 | return loc + scale*rk_gauss(state); |
|---|
| 105 | } |
|---|
| 106 | |
|---|
| 107 | double rk_standard_exponential(rk_state *state) |
|---|
| 108 | { |
|---|
| 109 | /* We use -log(1-U) since U is [0, 1) */ |
|---|
| 110 | return -log(1.0 - rk_double(state)); |
|---|
| 111 | } |
|---|
| 112 | |
|---|
| 113 | double rk_exponential(rk_state *state, double scale) |
|---|
| 114 | { |
|---|
| 115 | return scale * rk_standard_exponential(state); |
|---|
| 116 | } |
|---|
| 117 | |
|---|
| 118 | double rk_uniform(rk_state *state, double loc, double scale) |
|---|
| 119 | { |
|---|
| 120 | return loc + scale*rk_double(state); |
|---|
| 121 | } |
|---|
| 122 | |
|---|
| 123 | double rk_standard_gamma(rk_state *state, double shape) |
|---|
| 124 | { |
|---|
| 125 | double b, c; |
|---|
| 126 | double U, V, X, Y; |
|---|
| 127 | |
|---|
| 128 | if (shape == 1.0) |
|---|
| 129 | { |
|---|
| 130 | return rk_standard_exponential(state); |
|---|
| 131 | } |
|---|
| 132 | else if (shape < 1.0) |
|---|
| 133 | { |
|---|
| 134 | for (;;) |
|---|
| 135 | { |
|---|
| 136 | U = rk_double(state); |
|---|
| 137 | V = rk_standard_exponential(state); |
|---|
| 138 | if (U <= 1.0 - shape) |
|---|
| 139 | { |
|---|
| 140 | X = pow(U, 1./shape); |
|---|
| 141 | if (X <= V) |
|---|
| 142 | { |
|---|
| 143 | return X; |
|---|
| 144 | } |
|---|
| 145 | } |
|---|
| 146 | else |
|---|
| 147 | { |
|---|
| 148 | Y = -log((1-U)/shape); |
|---|
| 149 | X = pow(1.0 - shape + shape*Y, 1./shape); |
|---|
| 150 | if (X <= (V + Y)) |
|---|
| 151 | { |
|---|
| 152 | return X; |
|---|
| 153 | } |
|---|
| 154 | } |
|---|
| 155 | } |
|---|
| 156 | } |
|---|
| 157 | else |
|---|
| 158 | { |
|---|
| 159 | b = shape - 1./3.; |
|---|
| 160 | c = 1./sqrt(9*b); |
|---|
| 161 | for (;;) |
|---|
| 162 | { |
|---|
| 163 | do |
|---|
| 164 | { |
|---|
| 165 | X = rk_gauss(state); |
|---|
| 166 | V = 1.0 + c*X; |
|---|
| 167 | } while (V <= 0.0); |
|---|
| 168 | |
|---|
| 169 | V = V*V*V; |
|---|
| 170 | U = rk_double(state); |
|---|
| 171 | if (U < 1.0 - 0.0331*(X*X)*(X*X)) return (b*V); |
|---|
| 172 | if (log(U) < 0.5*X*X + b*(1. - V + log(V))) return (b*V); |
|---|
| 173 | } |
|---|
| 174 | } |
|---|
| 175 | } |
|---|
| 176 | |
|---|
| 177 | double rk_gamma(rk_state *state, double shape, double scale) |
|---|
| 178 | { |
|---|
| 179 | return scale * rk_standard_gamma(state, shape); |
|---|
| 180 | } |
|---|
| 181 | |
|---|
| 182 | double rk_beta(rk_state *state, double a, double b) |
|---|
| 183 | { |
|---|
| 184 | double Ga, Gb; |
|---|
| 185 | |
|---|
| 186 | if ((a <= 1.0) && (b <= 1.0)) |
|---|
| 187 | { |
|---|
| 188 | double U, V, X, Y; |
|---|
| 189 | /* Use Jonk's algorithm */ |
|---|
| 190 | |
|---|
| 191 | while (1) |
|---|
| 192 | { |
|---|
| 193 | U = rk_double(state); |
|---|
| 194 | V = rk_double(state); |
|---|
| 195 | X = pow(U, 1.0/a); |
|---|
| 196 | Y = pow(V, 1.0/b); |
|---|
| 197 | |
|---|
| 198 | if ((X + Y) <= 1.0) |
|---|
| 199 | { |
|---|
| 200 | return X / (X + Y); |
|---|
| 201 | } |
|---|
| 202 | } |
|---|
| 203 | } |
|---|
| 204 | else |
|---|
| 205 | { |
|---|
| 206 | Ga = rk_standard_gamma(state, a); |
|---|
| 207 | Gb = rk_standard_gamma(state, b); |
|---|
| 208 | return Ga/(Ga + Gb); |
|---|
| 209 | } |
|---|
| 210 | } |
|---|
| 211 | |
|---|
| 212 | double rk_chisquare(rk_state *state, double df) |
|---|
| 213 | { |
|---|
| 214 | return 2.0*rk_standard_gamma(state, df/2.0); |
|---|
| 215 | } |
|---|
| 216 | |
|---|
| 217 | double rk_noncentral_chisquare(rk_state *state, double df, double nonc) |
|---|
| 218 | { |
|---|
| 219 | double Chi2, N; |
|---|
| 220 | |
|---|
| 221 | Chi2 = rk_chisquare(state, df-1); |
|---|
| 222 | N = rk_gauss(state) + sqrt(nonc); |
|---|
| 223 | return Chi2 + N*N; |
|---|
| 224 | } |
|---|
| 225 | |
|---|
| 226 | double rk_f(rk_state *state, double dfnum, double dfden) |
|---|
| 227 | { |
|---|
| 228 | return ((rk_chisquare(state, dfnum) * dfden) / |
|---|
| 229 | (rk_chisquare(state, dfden) * dfnum)); |
|---|
| 230 | } |
|---|
| 231 | |
|---|
| 232 | double rk_noncentral_f(rk_state *state, double dfnum, double dfden, double nonc) |
|---|
| 233 | { |
|---|
| 234 | return ((rk_noncentral_chisquare(state, dfnum, nonc)*dfden) / |
|---|
| 235 | (rk_chisquare(state, dfden)*dfnum)); |
|---|
| 236 | } |
|---|
| 237 | |
|---|
| 238 | long rk_binomial_btpe(rk_state *state, long n, double p) |
|---|
| 239 | { |
|---|
| 240 | double r,q,fm,p1,xm,xl,xr,c,laml,lamr,p2,p3,p4; |
|---|
| 241 | double a,u,v,s,F,rho,t,A,nrq,x1,x2,f1,f2,z,z2,w,w2,x; |
|---|
| 242 | long m,y,k,i; |
|---|
| 243 | |
|---|
| 244 | if (!(state->has_binomial) || |
|---|
| 245 | (state->nsave != n) || |
|---|
| 246 | (state->psave != p)) |
|---|
| 247 | { |
|---|
| 248 | /* initialize */ |
|---|
| 249 | state->nsave = n; |
|---|
| 250 | state->psave = p; |
|---|
| 251 | state->has_binomial = 1; |
|---|
| 252 | state->r = r = min(p, 1.0-p); |
|---|
| 253 | state->q = q = 1.0 - r; |
|---|
| 254 | state->fm = fm = n*r+r; |
|---|
| 255 | state->m = m = (long)floor(state->fm); |
|---|
| 256 | state->p1 = p1 = floor(2.195*sqrt(n*r*q)-4.6*q) + 0.5; |
|---|
| 257 | state->xm = xm = m + 0.5; |
|---|
| 258 | state->xl = xl = xm - p1; |
|---|
| 259 | state->xr = xr = xm + p1; |
|---|
| 260 | state->c = c = 0.134 + 20.5/(15.3 + m); |
|---|
| 261 | a = (fm - xl)/(fm-xl*r); |
|---|
| 262 | state->laml = laml = a*(1.0 + a/2.0); |
|---|
| 263 | a = (xr - fm)/(xr*q); |
|---|
| 264 | state->lamr = lamr = a*(1.0 + a/2.0); |
|---|
| 265 | state->p2 = p2 = p1*(1.0 + 2.0*c); |
|---|
| 266 | state->p3 = p3 = p2 + c/laml; |
|---|
| 267 | state->p4 = p4 = p3 + c/lamr; |
|---|
| 268 | } |
|---|
| 269 | else |
|---|
| 270 | { |
|---|
| 271 | r = state->r; |
|---|
| 272 | q = state->q; |
|---|
| 273 | fm = state->fm; |
|---|
| 274 | m = state->m; |
|---|
| 275 | p1 = state->p1; |
|---|
| 276 | xm = state->xm; |
|---|
| 277 | xl = state->xl; |
|---|
| 278 | xr = state->xr; |
|---|
| 279 | c = state->c; |
|---|
| 280 | laml = state->laml; |
|---|
| 281 | lamr = state->lamr; |
|---|
| 282 | p2 = state->p2; |
|---|
| 283 | p3 = state->p3; |
|---|
| 284 | p4 = state->p4; |
|---|
| 285 | } |
|---|
| 286 | |
|---|
| 287 | /* sigh ... */ |
|---|
| 288 | Step10: |
|---|
| 289 | nrq = n*r*q; |
|---|
| 290 | u = rk_double(state)*p4; |
|---|
| 291 | v = rk_double(state); |
|---|
| 292 | if (u > p1) goto Step20; |
|---|
| 293 | y = (long)floor(xm - p1*v + u); |
|---|
| 294 | goto Step60; |
|---|
| 295 | |
|---|
| 296 | Step20: |
|---|
| 297 | if (u > p2) goto Step30; |
|---|
| 298 | x = xl + (u - p1)/c; |
|---|
| 299 | v = v*c + 1.0 - fabs(m - x + 0.5)/p1; |
|---|
| 300 | if (v > 1.0) goto Step10; |
|---|
| 301 | y = (long)floor(x); |
|---|
| 302 | goto Step50; |
|---|
| 303 | |
|---|
| 304 | Step30: |
|---|
| 305 | if (u > p3) goto Step40; |
|---|
| 306 | y = (long)floor(xl + log(v)/laml); |
|---|
| 307 | if (y < 0) goto Step10; |
|---|
| 308 | v = v*(u-p2)*laml; |
|---|
| 309 | goto Step50; |
|---|
| 310 | |
|---|
| 311 | Step40: |
|---|
| 312 | y = (int)floor(xr - log(v)/lamr); |
|---|
| 313 | if (y > n) goto Step10; |
|---|
| 314 | v = v*(u-p3)*lamr; |
|---|
| 315 | |
|---|
| 316 | Step50: |
|---|
| 317 | k = fabs(y - m); |
|---|
| 318 | if ((k > 20) && (k < ((nrq)/2.0 - 1))) goto Step52; |
|---|
| 319 | |
|---|
| 320 | s = r/q; |
|---|
| 321 | a = s*(n+1); |
|---|
| 322 | F = 1.0; |
|---|
| 323 | if (m < y) |
|---|
| 324 | { |
|---|
| 325 | for (i=m; i<=y; i++) |
|---|
| 326 | { |
|---|
| 327 | F *= (a/i - s); |
|---|
| 328 | } |
|---|
| 329 | } |
|---|
| 330 | else if (m > y) |
|---|
| 331 | { |
|---|
| 332 | for (i=y; i<=m; i++) |
|---|
| 333 | { |
|---|
| 334 | F /= (a/i - s); |
|---|
| 335 | } |
|---|
| 336 | } |
|---|
| 337 | else |
|---|
| 338 | { |
|---|
| 339 | if (v > F) goto Step10; |
|---|
| 340 | goto Step60; |
|---|
| 341 | } |
|---|
| 342 | |
|---|
| 343 | Step52: |
|---|
| 344 | rho = (k/(nrq))*((k*(k/3.0 + 0.625) + 0.16666666666666666)/nrq + 0.5); |
|---|
| 345 | t = -k*k/(2*nrq); |
|---|
| 346 | A = log(v); |
|---|
| 347 | if (A < (t - rho)) goto Step60; |
|---|
| 348 | if (A > (t + rho)) goto Step10; |
|---|
| 349 | |
|---|
| 350 | x1 = y+1; |
|---|
| 351 | f1 = m+1; |
|---|
| 352 | z = n+1-m; |
|---|
| 353 | w = n-y+1; |
|---|
| 354 | x2 = x1*x1; |
|---|
| 355 | f2 = f1*f1; |
|---|
| 356 | z2 = z*z; |
|---|
| 357 | w2 = w*w; |
|---|
| 358 | if (A > (xm*log(f1/x1) |
|---|
| 359 | + (n-m+0.5)*log(z/w) |
|---|
| 360 | + (y-m)*log(w*r/(x1*q)) |
|---|
| 361 | + (13680.-(462.-(132.-(99.-140./f2)/f2)/f2)/f2)/f1/166320. |
|---|
| 362 | + (13680.-(462.-(132.-(99.-140./z2)/z2)/z2)/z2)/z/166320. |
|---|
| 363 | + (13680.-(462.-(132.-(99.-140./x2)/x2)/x2)/x2)/x1/166320. |
|---|
| 364 | + (13680.-(462.-(132.-(99.-140./w2)/w2)/w2)/w2)/w/166320.)) |
|---|
| 365 | { |
|---|
| 366 | goto Step10; |
|---|
| 367 | } |
|---|
| 368 | |
|---|
| 369 | Step60: |
|---|
| 370 | if (p > 0.5) |
|---|
| 371 | { |
|---|
| 372 | y = n - y; |
|---|
| 373 | } |
|---|
| 374 | |
|---|
| 375 | return y; |
|---|
| 376 | } |
|---|
| 377 | |
|---|
| 378 | long rk_binomial_inversion(rk_state *state, long n, double p) |
|---|
| 379 | { |
|---|
| 380 | double q, qn, np, px, U; |
|---|
| 381 | long X, bound; |
|---|
| 382 | |
|---|
| 383 | if (!(state->has_binomial) || |
|---|
| 384 | (state->nsave != n) || |
|---|
| 385 | (state->psave != p)) |
|---|
| 386 | { |
|---|
| 387 | state->nsave = n; |
|---|
| 388 | state->psave = p; |
|---|
| 389 | state->has_binomial = 1; |
|---|
| 390 | state->q = q = 1.0 - p; |
|---|
| 391 | state->r = qn = exp(n * log(q)); |
|---|
| 392 | state->c = np = n*p; |
|---|
| 393 | state->m = bound = min(n, np + 10.0*sqrt(np*q + 1)); |
|---|
| 394 | } else |
|---|
| 395 | { |
|---|
| 396 | q = state->q; |
|---|
| 397 | qn = state->r; |
|---|
| 398 | np = state->c; |
|---|
| 399 | bound = state->m; |
|---|
| 400 | } |
|---|
| 401 | X = 0; |
|---|
| 402 | px = qn; |
|---|
| 403 | U = rk_double(state); |
|---|
| 404 | while (U > px) |
|---|
| 405 | { |
|---|
| 406 | X++; |
|---|
| 407 | if (X > bound) |
|---|
| 408 | { |
|---|
| 409 | X = 0; |
|---|
| 410 | px = qn; |
|---|
| 411 | U = rk_double(state); |
|---|
| 412 | } else |
|---|
| 413 | { |
|---|
| 414 | U -= px; |
|---|
| 415 | px = ((n-X+1) * p * px)/(X*q); |
|---|
| 416 | } |
|---|
| 417 | } |
|---|
| 418 | return X; |
|---|
| 419 | } |
|---|
| 420 | |
|---|
| 421 | long rk_binomial(rk_state *state, long n, double p) |
|---|
| 422 | { |
|---|
| 423 | double q; |
|---|
| 424 | |
|---|
| 425 | if (p <= 0.5) |
|---|
| 426 | { |
|---|
| 427 | if (p*n <= 30.0) |
|---|
| 428 | { |
|---|
| 429 | return rk_binomial_inversion(state, n, p); |
|---|
| 430 | } |
|---|
| 431 | else |
|---|
| 432 | { |
|---|
| 433 | return rk_binomial_btpe(state, n, p); |
|---|
| 434 | } |
|---|
| 435 | } |
|---|
| 436 | else |
|---|
| 437 | { |
|---|
| 438 | q = 1.0-p; |
|---|
| 439 | if (q*n <= 30.0) |
|---|
| 440 | { |
|---|
| 441 | return n - rk_binomial_inversion(state, n, q); |
|---|
| 442 | } |
|---|
| 443 | else |
|---|
| 444 | { |
|---|
| 445 | return n - rk_binomial_btpe(state, n, q); |
|---|
| 446 | } |
|---|
| 447 | } |
|---|
| 448 | |
|---|
| 449 | } |
|---|
| 450 | |
|---|
| 451 | long rk_negative_binomial(rk_state *state, double n, double p) |
|---|
| 452 | { |
|---|
| 453 | double Y; |
|---|
| 454 | |
|---|
| 455 | Y = rk_gamma(state, n, (1-p)/p); |
|---|
| 456 | return rk_poisson(state, Y); |
|---|
| 457 | } |
|---|
| 458 | |
|---|
| 459 | long rk_poisson_mult(rk_state *state, double lam) |
|---|
| 460 | { |
|---|
| 461 | long X; |
|---|
| 462 | double prod, U, enlam; |
|---|
| 463 | |
|---|
| 464 | enlam = exp(-lam); |
|---|
| 465 | X = 0; |
|---|
| 466 | prod = 1.0; |
|---|
| 467 | while (1) |
|---|
| 468 | { |
|---|
| 469 | U = rk_double(state); |
|---|
| 470 | prod *= U; |
|---|
| 471 | if (prod > enlam) |
|---|
| 472 | { |
|---|
| 473 | X += 1; |
|---|
| 474 | } |
|---|
| 475 | else |
|---|
| 476 | { |
|---|
| 477 | return X; |
|---|
| 478 | } |
|---|
| 479 | } |
|---|
| 480 | } |
|---|
| 481 | |
|---|
| 482 | #define LS2PI 0.91893853320467267 |
|---|
| 483 | #define TWELFTH 0.083333333333333333333333 |
|---|
| 484 | long rk_poisson_ptrs(rk_state *state, double lam) |
|---|
| 485 | { |
|---|
| 486 | long k; |
|---|
| 487 | double U, V, slam, loglam, a, b, invalpha, vr, us; |
|---|
| 488 | |
|---|
| 489 | slam = sqrt(lam); |
|---|
| 490 | loglam = log(lam); |
|---|
| 491 | b = 0.931 + 2.53*slam; |
|---|
| 492 | a = -0.059 + 0.02483*b; |
|---|
| 493 | invalpha = 1.1239 + 1.1328/(b-3.4); |
|---|
| 494 | vr = 0.9277 - 3.6224/(b-2); |
|---|
| 495 | |
|---|
| 496 | while (1) |
|---|
| 497 | { |
|---|
| 498 | U = rk_double(state) - 0.5; |
|---|
| 499 | V = rk_double(state); |
|---|
| 500 | us = 0.5 - fabs(U); |
|---|
| 501 | k = (long)floor((2*a/us + b)*U + lam + 0.43); |
|---|
| 502 | if ((us >= 0.07) && (V <= vr)) |
|---|
| 503 | { |
|---|
| 504 | return k; |
|---|
| 505 | } |
|---|
| 506 | if ((k < 0) || |
|---|
| 507 | ((us < 0.013) && (V > us))) |
|---|
| 508 | { |
|---|
| 509 | continue; |
|---|
| 510 | } |
|---|
| 511 | if ((log(V) + log(invalpha) - log(a/(us*us)+b)) <= |
|---|
| 512 | (-lam + k*loglam - loggam(k+1))) |
|---|
| 513 | { |
|---|
| 514 | return k; |
|---|
| 515 | } |
|---|
| 516 | |
|---|
| 517 | |
|---|
| 518 | } |
|---|
| 519 | |
|---|
| 520 | } |
|---|
| 521 | |
|---|
| 522 | long rk_poisson(rk_state *state, double lam) |
|---|
| 523 | { |
|---|
| 524 | if (lam >= 10) |
|---|
| 525 | { |
|---|
| 526 | return rk_poisson_ptrs(state, lam); |
|---|
| 527 | } |
|---|
| 528 | else if (lam == 0) |
|---|
| 529 | { |
|---|
| 530 | return 0; |
|---|
| 531 | } |
|---|
| 532 | else |
|---|
| 533 | { |
|---|
| 534 | return rk_poisson_mult(state, lam); |
|---|
| 535 | } |
|---|
| 536 | } |
|---|
| 537 | |
|---|
| 538 | double rk_standard_cauchy(rk_state *state) |
|---|
| 539 | { |
|---|
| 540 | return rk_gauss(state) / rk_gauss(state); |
|---|
| 541 | } |
|---|
| 542 | |
|---|
| 543 | double rk_standard_t(rk_state *state, double df) |
|---|
| 544 | { |
|---|
| 545 | double N, G, X; |
|---|
| 546 | |
|---|
| 547 | N = rk_gauss(state); |
|---|
| 548 | G = rk_standard_gamma(state, df/2); |
|---|
| 549 | X = sqrt(df/2)*N/sqrt(G); |
|---|
| 550 | return X; |
|---|
| 551 | } |
|---|
| 552 | |
|---|
| 553 | /* Uses the rejection algorithm compared against the wrapped Cauchy |
|---|
| 554 | distribution suggested by Best and Fisher and documented in |
|---|
| 555 | Chapter 9 of Luc's Non-Uniform Random Variate Generation. |
|---|
| 556 | http://cg.scs.carleton.ca/~luc/rnbookindex.html |
|---|
| 557 | (but corrected to match the algorithm in R and Python) |
|---|
| 558 | */ |
|---|
| 559 | double rk_vonmises(rk_state *state, double mu, double kappa) |
|---|
| 560 | { |
|---|
| 561 | double r, rho, s; |
|---|
| 562 | double U, V, W, Y, Z; |
|---|
| 563 | double result, mod; |
|---|
| 564 | int neg; |
|---|
| 565 | |
|---|
| 566 | if (kappa < 1e-8) |
|---|
| 567 | { |
|---|
| 568 | return M_PI * (2*rk_double(state)-1); |
|---|
| 569 | } |
|---|
| 570 | else |
|---|
| 571 | { |
|---|
| 572 | r = 1 + sqrt(1 + 4*kappa*kappa); |
|---|
| 573 | rho = (r - sqrt(2*r))/(2*kappa); |
|---|
| 574 | s = (1 + rho*rho)/(2*rho); |
|---|
| 575 | |
|---|
| 576 | while (1) |
|---|
| 577 | { |
|---|
| 578 | U = rk_double(state); |
|---|
| 579 | Z = cos(M_PI*U); |
|---|
| 580 | W = (1 + s*Z)/(s + Z); |
|---|
| 581 | Y = kappa * (s - W); |
|---|
| 582 | V = rk_double(state); |
|---|
| 583 | if ((Y*(2-Y) - V >= 0) || (log(Y/V)+1 - Y >= 0)) |
|---|
| 584 | { |
|---|
| 585 | break; |
|---|
| 586 | } |
|---|
| 587 | } |
|---|
| 588 | |
|---|
| 589 | U = rk_double(state); |
|---|
| 590 | |
|---|
| 591 | result = acos(W); |
|---|
| 592 | if (U < 0.5) |
|---|
| 593 | { |
|---|
| 594 | result = -result; |
|---|
| 595 | } |
|---|
| 596 | result += mu; |
|---|
| 597 | neg = (result < 0); |
|---|
| 598 | mod = fabs(result); |
|---|
| 599 | mod = (fmod(mod+M_PI, 2*M_PI)-M_PI); |
|---|
| 600 | if (neg) |
|---|
| 601 | { |
|---|
| 602 | mod *= -1; |
|---|
| 603 | } |
|---|
| 604 | |
|---|
| 605 | return mod; |
|---|
| 606 | } |
|---|
| 607 | } |
|---|
| 608 | |
|---|
| 609 | double rk_pareto(rk_state *state, double a) |
|---|
| 610 | { |
|---|
| 611 | return exp(rk_standard_exponential(state)/a) - 1; |
|---|
| 612 | } |
|---|
| 613 | |
|---|
| 614 | double rk_weibull(rk_state *state, double a) |
|---|
| 615 | { |
|---|
| 616 | return pow(rk_standard_exponential(state), 1./a); |
|---|
| 617 | } |
|---|
| 618 | |
|---|
| 619 | double rk_power(rk_state *state, double a) |
|---|
| 620 | { |
|---|
| 621 | return pow(1 - exp(-rk_standard_exponential(state)), 1./a); |
|---|
| 622 | } |
|---|
| 623 | |
|---|
| 624 | double rk_laplace(rk_state *state, double loc, double scale) |
|---|
| 625 | { |
|---|
| 626 | double U; |
|---|
| 627 | |
|---|
| 628 | U = rk_double(state); |
|---|
| 629 | if (U < 0.5) |
|---|
| 630 | { |
|---|
| 631 | U = loc + scale * log(U + U); |
|---|
| 632 | } else |
|---|
| 633 | { |
|---|
| 634 | U = loc - scale * log(2.0 - U - U); |
|---|
| 635 | } |
|---|
| 636 | return U; |
|---|
| 637 | } |
|---|
| 638 | |
|---|
| 639 | double rk_gumbel(rk_state *state, double loc, double scale) |
|---|
| 640 | { |
|---|
| 641 | double U; |
|---|
| 642 | |
|---|
| 643 | U = 1.0 - rk_double(state); |
|---|
| 644 | return loc - scale * log(-log(U)); |
|---|
| 645 | } |
|---|
| 646 | |
|---|
| 647 | double rk_logistic(rk_state *state, double loc, double scale) |
|---|
| 648 | { |
|---|
| 649 | double U; |
|---|
| 650 | |
|---|
| 651 | U = rk_double(state); |
|---|
| 652 | return loc + scale * log(U/(1.0 - U)); |
|---|
| 653 | } |
|---|
| 654 | |
|---|
| 655 | double rk_lognormal(rk_state *state, double mean, double sigma) |
|---|
| 656 | { |
|---|
| 657 | return exp(rk_normal(state, mean, sigma)); |
|---|
| 658 | } |
|---|
| 659 | |
|---|
| 660 | double rk_rayleigh(rk_state *state, double mode) |
|---|
| 661 | { |
|---|
| 662 | return mode*sqrt(-2.0 * log(1.0 - rk_double(state))); |
|---|
| 663 | } |
|---|
| 664 | |
|---|
| 665 | double rk_wald(rk_state *state, double mean, double scale) |
|---|
| 666 | { |
|---|
| 667 | double U, X, Y; |
|---|
| 668 | double mu_2l; |
|---|
| 669 | |
|---|
| 670 | mu_2l = mean / (2*scale); |
|---|
| 671 | Y = rk_gauss(state); |
|---|
| 672 | Y = mean*Y*Y; |
|---|
| 673 | X = mean + mu_2l*(Y - sqrt(4*scale*Y + Y*Y)); |
|---|
| 674 | U = rk_double(state); |
|---|
| 675 | if (U <= mean/(mean+X)) |
|---|
| 676 | { |
|---|
| 677 | return X; |
|---|
| 678 | } else |
|---|
| 679 | { |
|---|
| 680 | return mean*mean/X; |
|---|
| 681 | } |
|---|
| 682 | } |
|---|
| 683 | |
|---|
| 684 | long rk_zipf(rk_state *state, double a) |
|---|
| 685 | { |
|---|
| 686 | double T, U, V; |
|---|
| 687 | long X; |
|---|
| 688 | double am1, b; |
|---|
| 689 | |
|---|
| 690 | am1 = a - 1.0; |
|---|
| 691 | b = pow(2.0, am1); |
|---|
| 692 | do |
|---|
| 693 | { |
|---|
| 694 | U = 1.0-rk_double(state); |
|---|
| 695 | V = rk_double(state); |
|---|
| 696 | X = (long)floor(pow(U, -1.0/am1)); |
|---|
| 697 | /* The real result may be above what can be represented in a signed |
|---|
| 698 | * long. It will get casted to -sys.maxint-1. Since this is |
|---|
| 699 | * a straightforward rejection algorithm, we can just reject this value |
|---|
| 700 | * in the rejection condition below. This function then models a Zipf |
|---|
| 701 | * distribution truncated to sys.maxint. |
|---|
| 702 | */ |
|---|
| 703 | T = pow(1.0 + 1.0/X, am1); |
|---|
| 704 | } while (((V*X*(T-1.0)/(b-1.0)) > (T/b)) || X < 1); |
|---|
| 705 | return X; |
|---|
| 706 | } |
|---|
| 707 | |
|---|
| 708 | long rk_geometric_search(rk_state *state, double p) |
|---|
| 709 | { |
|---|
| 710 | double U; |
|---|
| 711 | long X; |
|---|
| 712 | double sum, prod, q; |
|---|
| 713 | |
|---|
| 714 | X = 1; |
|---|
| 715 | sum = prod = p; |
|---|
| 716 | q = 1.0 - p; |
|---|
| 717 | U = rk_double(state); |
|---|
| 718 | while (U > sum) |
|---|
| 719 | { |
|---|
| 720 | prod *= q; |
|---|
| 721 | sum += prod; |
|---|
| 722 | X++; |
|---|
| 723 | } |
|---|
| 724 | return X; |
|---|
| 725 | } |
|---|
| 726 | |
|---|
| 727 | long rk_geometric_inversion(rk_state *state, double p) |
|---|
| 728 | { |
|---|
| 729 | return (long)ceil(log(1.0-rk_double(state))/log(1.0-p)); |
|---|
| 730 | } |
|---|
| 731 | |
|---|
| 732 | long rk_geometric(rk_state *state, double p) |
|---|
| 733 | { |
|---|
| 734 | if (p >= 0.333333333333333333333333) |
|---|
| 735 | { |
|---|
| 736 | return rk_geometric_search(state, p); |
|---|
| 737 | } else |
|---|
| 738 | { |
|---|
| 739 | return rk_geometric_inversion(state, p); |
|---|
| 740 | } |
|---|
| 741 | } |
|---|
| 742 | |
|---|
| 743 | long rk_hypergeometric_hyp(rk_state *state, long good, long bad, long sample) |
|---|
| 744 | { |
|---|
| 745 | long d1, K, Z; |
|---|
| 746 | double d2, U, Y; |
|---|
| 747 | |
|---|
| 748 | d1 = bad + good - sample; |
|---|
| 749 | d2 = (double)min(bad, good); |
|---|
| 750 | |
|---|
| 751 | Y = d2; |
|---|
| 752 | K = sample; |
|---|
| 753 | while (Y > 0.0) |
|---|
| 754 | { |
|---|
| 755 | U = rk_double(state); |
|---|
| 756 | Y -= (long)floor(U + Y/(d1 + K)); |
|---|
| 757 | K--; |
|---|
| 758 | if (K == 0) break; |
|---|
| 759 | } |
|---|
| 760 | Z = (long)(d2 - Y); |
|---|
| 761 | if (good > bad) Z = sample - Z; |
|---|
| 762 | return Z; |
|---|
| 763 | } |
|---|
| 764 | |
|---|
| 765 | /* D1 = 2*sqrt(2/e) */ |
|---|
| 766 | /* D2 = 3 - 2*sqrt(3/e) */ |
|---|
| 767 | #define D1 1.7155277699214135 |
|---|
| 768 | #define D2 0.8989161620588988 |
|---|
| 769 | long rk_hypergeometric_hrua(rk_state *state, long good, long bad, long sample) |
|---|
| 770 | { |
|---|
| 771 | long mingoodbad, maxgoodbad, popsize, m, d9; |
|---|
| 772 | double d4, d5, d6, d7, d8, d10, d11; |
|---|
| 773 | long Z; |
|---|
| 774 | double T, W, X, Y; |
|---|
| 775 | |
|---|
| 776 | mingoodbad = min(good, bad); |
|---|
| 777 | popsize = good + bad; |
|---|
| 778 | maxgoodbad = max(good, bad); |
|---|
| 779 | m = min(sample, popsize - sample); |
|---|
| 780 | d4 = ((double)mingoodbad) / popsize; |
|---|
| 781 | d5 = 1.0 - d4; |
|---|
| 782 | d6 = m*d4 + 0.5; |
|---|
| 783 | d7 = sqrt((popsize - m) * sample * d4 *d5 / (popsize-1) + 0.5); |
|---|
| 784 | d8 = D1*d7 + D2; |
|---|
| 785 | d9 = (long)floor((double)((m+1)*(mingoodbad+1))/(popsize+2)); |
|---|
| 786 | d10 = (loggam(d9+1) + loggam(mingoodbad-d9+1) + loggam(m-d9+1) + |
|---|
| 787 | loggam(maxgoodbad-m+d9+1)); |
|---|
| 788 | d11 = min(min(m, mingoodbad)+1.0, floor(d6+16*d7)); |
|---|
| 789 | /* 16 for 16-decimal-digit precision in D1 and D2 */ |
|---|
| 790 | |
|---|
| 791 | while (1) |
|---|
| 792 | { |
|---|
| 793 | X = rk_double(state); |
|---|
| 794 | Y = rk_double(state); |
|---|
| 795 | W = d6 + d8*(Y- 0.5)/X; |
|---|
| 796 | |
|---|
| 797 | /* fast rejection: */ |
|---|
| 798 | if ((W < 0.0) || (W >= d11)) continue; |
|---|
| 799 | |
|---|
| 800 | Z = (long)floor(W); |
|---|
| 801 | T = d10 - (loggam(Z+1) + loggam(mingoodbad-Z+1) + loggam(m-Z+1) + |
|---|
| 802 | loggam(maxgoodbad-m+Z+1)); |
|---|
| 803 | |
|---|
| 804 | /* fast acceptance: */ |
|---|
| 805 | if ((X*(4.0-X)-3.0) <= T) break; |
|---|
| 806 | |
|---|
| 807 | /* fast rejection: */ |
|---|
| 808 | if (X*(X-T) >= 1) continue; |
|---|
| 809 | |
|---|
| 810 | if (2.0*log(X) <= T) break; /* acceptance */ |
|---|
| 811 | } |
|---|
| 812 | |
|---|
| 813 | /* this is a correction to HRUA* by Ivan Frohne in rv.py */ |
|---|
| 814 | if (good > bad) Z = m - Z; |
|---|
| 815 | |
|---|
| 816 | /* another fix from rv.py to allow sample to exceed popsize/2 */ |
|---|
| 817 | if (m < sample) Z = good - Z; |
|---|
| 818 | |
|---|
| 819 | return Z; |
|---|
| 820 | } |
|---|
| 821 | #undef D1 |
|---|
| 822 | #undef D2 |
|---|
| 823 | |
|---|
| 824 | long rk_hypergeometric(rk_state *state, long good, long bad, long sample) |
|---|
| 825 | { |
|---|
| 826 | if (sample > 10) |
|---|
| 827 | { |
|---|
| 828 | return rk_hypergeometric_hrua(state, good, bad, sample); |
|---|
| 829 | } else |
|---|
| 830 | { |
|---|
| 831 | return rk_hypergeometric_hyp(state, good, bad, sample); |
|---|
| 832 | } |
|---|
| 833 | } |
|---|
| 834 | |
|---|
| 835 | double rk_triangular(rk_state *state, double left, double mode, double right) |
|---|
| 836 | { |
|---|
| 837 | double base, leftbase, ratio, leftprod, rightprod; |
|---|
| 838 | double U; |
|---|
| 839 | |
|---|
| 840 | base = right - left; |
|---|
| 841 | leftbase = mode - left; |
|---|
| 842 | ratio = leftbase / base; |
|---|
| 843 | leftprod = leftbase*base; |
|---|
| 844 | rightprod = (right - mode)*base; |
|---|
| 845 | |
|---|
| 846 | U = rk_double(state); |
|---|
| 847 | if (U <= ratio) |
|---|
| 848 | { |
|---|
| 849 | return left + sqrt(U*leftprod); |
|---|
| 850 | } else |
|---|
| 851 | { |
|---|
| 852 | return right - sqrt((1.0 - U) * rightprod); |
|---|
| 853 | } |
|---|
| 854 | } |
|---|
| 855 | |
|---|
| 856 | long rk_logseries(rk_state *state, double p) |
|---|
| 857 | { |
|---|
| 858 | double q, r, U, V; |
|---|
| 859 | long result; |
|---|
| 860 | |
|---|
| 861 | r = log(1.0 - p); |
|---|
| 862 | |
|---|
| 863 | while (1) { |
|---|
| 864 | V = rk_double(state); |
|---|
| 865 | if (V >= p) { |
|---|
| 866 | return 1; |
|---|
| 867 | } |
|---|
| 868 | U = rk_double(state); |
|---|
| 869 | q = 1.0 - exp(r*U); |
|---|
| 870 | if (V <= q*q) { |
|---|
| 871 | result = (long)floor(1 + log(V)/log(q)); |
|---|
| 872 | if (result < 1) { |
|---|
| 873 | continue; |
|---|
| 874 | } |
|---|
| 875 | else { |
|---|
| 876 | return result; |
|---|
| 877 | } |
|---|
| 878 | } |
|---|
| 879 | if (V >= q) { |
|---|
| 880 | return 1; |
|---|
| 881 | } |
|---|
| 882 | return 2; |
|---|
| 883 | } |
|---|
| 884 | } |
|---|