#### d|ensity #### p|robability (cumulative) #### q|uantile #### r|andom number generation #### #### Functions for ``d/p/q/r'' === "regression" tests (no *.Rout.save) options(warn = 2) ## ======== No warnings, unless explicitly asserted via assertWarning <- tools::assertWarning as.nan <- function(x) { x[is.na(x) & !is.nan(x)] <- NaN ; x } ###-- these are identical in ./arith-true.R ["fixme": use source(..)] ## opt.conformance <- 0 (sysinf <- Sys.info()) Lnx <- sysinf[["sysname"]] == "Linux" isMac <- sysinf[["sysname"]] == "Darwin" arch <- sysinf[["machine"]] onWindows <- .Platform$OS.type == "windows" b64 <- .Machine$sizeof.pointer >= 8 # 64 (or more) bits str(.Machine[grep("^sizeof", names(.Machine))]) ## also differentiate long-double.. (usingMKL <- grepl("/(lib)?mkl", La_library(), ignore.case=TRUE)) x86 <- arch == "x86_64" options(rErr.eps = 1e-30) rErr <- function(approx, true, eps = getOption("rErr.eps", 1e-30)) { ifelse(Mod(true) >= eps, 1 - approx / true, # relative error true - approx) # absolute error (e.g. when true=0) } ## Numerical equality: Here want "rel.error" almost always: All.eq <- function(x,y) { all.equal.numeric(x,y, tolerance = 64*.Machine$double.eps, scale = max(0, mean(abs(x), na.rm=TRUE))) } if(!interactive()) set.seed(123) .ptime <- proc.time() ## The prefixes of ALL the PDQ & R functions PDQRinteg <- c("binom", "geom", "hyper", "nbinom", "pois","signrank","wilcox") PDQR <- c(PDQRinteg, "beta", "cauchy", "chisq", "exp", "f", "gamma", "lnorm", "logis", "norm", "t","unif","weibull") PQonly <- c("tukey") ### (Extreme) tail tests added more recently: All.eq(1, -1e-17/ pexp(qexp(-1e-17, log=TRUE),log=TRUE)) abs(pgamma(30,100, lower=FALSE, log=TRUE) + 7.3384686328784e-24) < 1e-36 All.eq(1, pcauchy(-1e20) / 3.18309886183791e-21) All.eq(1, pcauchy(+1e15, log=TRUE) / -3.18309886183791e-16)## PR#6756 x <- 10^(ex <- c(1,2,5*(1:5),50,100,200,300,Inf)) for(a in x[ex > 10]) ## improve pt() : cbind(x,t= pt(-x, df=1), C=pcauchy(-x)) stopifnot(all.equal(pt(-a, df=1), pcauchy(-a), tolerance = 1e-15)) ## for PR#7902: ex <- -c(rev(1/x), ex) All.eq(-x, qcauchy(pcauchy(-x))) All.eq(+x, qcauchy(pcauchy(+x, log=TRUE), log=TRUE)) All.eq(1/x, pcauchy(qcauchy(1/x))) All.eq(ex, pcauchy(qcauchy(ex, log=TRUE), log=TRUE)) II <- c(-Inf,Inf) stopifnot(pcauchy(II) == 0:1, qcauchy(0:1) == II, pcauchy(II, log=TRUE) == c(-Inf,0), qcauchy(c(-Inf,0), log=TRUE) == II) ## PR#15521 : p <- 1 - 1/4096 stopifnot(all.equal(qcauchy(p), 1303.7970381453319163, tolerance = 1e-14)) pr <- 1e-23 ## PR#6757 stopifnot(all.equal(pr^ 12, pbinom(11, 12, prob= pr,lower=FALSE), tolerance = 1e-12, scale= 1e-270)) ## pbinom(.) gave 0 in R 1.9.0 pp <- 1e-17 ## PR#6792 stopifnot(all.equal(2*pp, pgeom(1, pp), scale= 1e-20)) ## pgeom(.) gave 0 in R 1.9.0 x <- 10^(100:295) sapply(c(1e-250, 1e-25, 0.9, 1.1, 101, 1e10, 1e100), function(shape) All.eq(-x, pgamma(x, shape=shape, lower=FALSE, log=TRUE))) x <- 2^(-1022:-900) ## where all completely off in R 2.0.1 all.equal(pgamma(x, 10, log = TRUE) - 10*log(x), rep(-15.104412573076, length(x)), tolerance = 1e-12)# 3.984e-14 (i386) all.equal(pgamma(x, 0.1, log = TRUE) - 0.1*log(x), rep(0.0498724412598364, length(x)), tolerance = 1e-13)# 7e-16 (i386) All.eq(dpois( 10*1:2, 3e-308, log=TRUE), c(-7096.08037610806, -14204.2875435307)) All.eq(dpois(1e20, 1e-290, log=TRUE), -7.12801378828154e+22) ## all gave -Inf in R 2.0.1 x <- c(outer(1:12, 10^c(-3:2, 6:9, 10*(2:30)))) for(nu in c(.75, 1.2, 4.5, 999, 1e50)) { lfx <- dt(x, df=nu, log=TRUE) stopifnot(is.finite(lfx), All.eq(exp(lfx), dt(x, df=nu))) }## dt(1e160, 1.2, log=TRUE) was -Inf up to R 2.15.2 ## pf() with large df1 or df2 ## (was said to be PR#7099, but that is about non-central pchisq) nu <- 2^seq(25, 34, 0.5) target <- pchisq(1, 1) # 0.682... y <- pf(1, 1, nu) stopifnot(All.eq(pf(1, 1, Inf), target), diff(c(y, target)) > 0, # i.e. pf(1, 1, *) is monotone increasing abs(y[1] - (target - 7.21129e-9)) < 1e-11) # computed value ## non-monotone in R <= 2.1.0 stopifnot(pgamma(Inf, 1.1) == 1) ## didn't not terminate in R 2.1.x (only) ## qgamma(q, *) should give {0,Inf} for q={0,1} sh <- c(1.1, 0.5, 0.2, 0.15, 1e-2, 1e-10) stopifnot(Inf == qgamma(1, sh)) stopifnot(0 == qgamma(0, sh)) ## the first gave Inf, NaN, and 99.425 in R 2.1.1 and earlier ## In extreme left tail {PR#11030} p <- 10:123*1e-12 qg <- qgamma(p, shape=19) qg2<- qgamma(1:100 * 1e-9, shape=11) stopifnot(diff(qg, diff=2) < -6e-6, diff(qg2,diff=2) < -6e-6, abs(1 - pgamma(qg, 19)/ p) < 1e-13, All.eq(qg [1], 2.35047385139143), All.eq(qg2[30], 1.11512318734547)) ## was non-continuous in R 2.6.2 and earlier f2 <- c(0.5, 1:4) stopifnot(df(0, 1, f2) == Inf, df(0, 2, f2) == 1, df(0, 3, f2) == 0) ## only the last one was ok in R 2.2.1 and earlier x0 <- -2 * 10^-c(22,10,7,5) # ==> d*() warns about non-integer: assertWarning(fx0 <- dbinom(x0, size = 3, prob = 0.1)) stopifnot(fx0 == 0, pbinom(x0, size = 3, prob = 0.1) == 0) ## very small negatives were rounded to 0 in R 2.2.1 and earlier ## dbeta(*, ncp): db.x <- c(0, 5, 80, 405, 1280, 3125, 6480, 12005, 20480, 32805, 50000, 73205, 103680, 142805, 192080, 253125, 327680) a <- rlnorm(100) stopifnot(All.eq(a, dbeta(0, 1, a, ncp=0)), dbeta(0, 0.9, 2.2, ncp = c(0, a)) == Inf, All.eq(65536 * dbeta(0:16/16, 5,1), db.x), All.eq(exp(16 * log(2) + dbeta(0:16/16, 5,1, log=TRUE)), db.x) ) ## the first gave 0, the 2nd NaN in R <= 2.3.0; others use 'TRUE' values stopifnot(all.equal(dbeta(0.8, 0.5, 5, ncp=1000),# was way too small in R <= 2.6.2 3.001852308909e-35), all.equal(1, integrate(dbeta, 0,1, 0.8, 0.5, ncp=1000)$value, tolerance = 1e-4), all.equal(1, integrate(dbeta, 0,1, 0.5, 200, ncp=720)$value), all.equal(1, integrate(dbeta, 0,1, 125, 200, ncp=2000)$value) ) ## df(*, ncp): x <- seq(0, 10, length=101) h <- 1e-7 dx.h <- (pf(x+h, 7, 5, ncp= 2.5) - pf(x-h, 7, 5, ncp= 2.5)) / (2*h) stopifnot(all.equal(dx.h, df(x, 7, 5, ncp= 2.5), tolerance = 1e-6),# (1.50 | 1.65)e-8 All.eq(df(0, 2, 4, ncp=x), df(1e-300, 2, 4, ncp=x)) ) ## qt(p ~ 0, df=1) - PR#9804 p <- 10^(-10:-20) qtp <- qt(p, df = 1) ## relative error: table(rerr <- abs(1 - p / pt(qtp, df=1))) # 64-bit F34 : 1 x 0, 10 x 2.22e-16 stopifnot(rerr < 1e-14) ## Similarly for df = 2 --- both for p ~ 0 *and* p ~ 1/2 ## P ~ 0 stopifnot(all.equal(qt(-740, df=2, log=TRUE), -exp(370)/sqrt(2))) ## P ~ 1 (=> p ~ 0.5): p.5 <- 0.5 + 2^(-5*(5:8)) stopifnot(all.equal(qt(p.5, df = 2), c(8.429369702179e-08, 2.634178031931e-09, 8.231806349784e-11, 2.572439484308e-12))) ## qt(, log = TRUE) is now more finite and monotone (again!): stopifnot(all.equal(qt(-1000, df = 4, log=TRUE), -4.930611e108, tolerance = 1e-6)) qtp <- qt(-(20:850), df=1.2, log=TRUE, lower=FALSE) ##almost: stopifnot(all(abs(5/6 - diff(log(qtp))) < 1e-11)) stopifnot(abs(5/6 - quantile(diff(log(qtp)), pr=c(0,0.995))) < 1e-11) ## close to df=1 (where Taylor steps are important!): stopifnot(all.equal(-20, pt(qt(-20, df=1.02, log=TRUE), df=1.02, log=TRUE), tolerance = 1e-12), diff(lq <- log(qt(-2^-(10:600), df=1.1, log=TRUE))) > 0.6) lq1 <- log(qt(-2^-(20:600), df=1, log=TRUE)) lq2 <- log(qt(-2^-(20:600), df=2, log=TRUE)) stopifnot(mean(abs(diff(lq1) - log(2) )) < 1e-8, mean(abs(diff(lq2) - log(sqrt(2)))) < 4e-8) ## Case, where log.p=TRUE was fine, but log.p=FALSE (default) gave NaN: lp <- 40:406 stopifnot(all.equal(lp, -pt(qt(exp(-lp), 1.2), 1.2, log=TRUE), tolerance = 4e-16)) ## Other log.p cases, gave NaN (all but 1.1) in R <= 4.2.1, PR#18360 [NB: *still* inaccurate: tol=0.2] q <- exp(seq(200, 500, by=1/2)) for(df in c(1.001, 1 + (1:10)/100)) { pq <- pt(q, df = df, log = TRUE) qpq <- qt(pq, df = df, log = TRUE) cat("df = ", df, ": all.equal(., tol=0): "); print(all.equal(q,qpq, tol=0)) # ~0.17! ## plot(lp, 1-qpq/q, main=paste0("relErr(qt(., df=",df,"))"), type="l") stopifnot(all.equal(q,qpq, tol = 0.2)) # Lnx 64b: 1.001 shows 0.179 if(any(ina <- is.na(qpq))) { cat("NaN in q-range: [", range(q[ina]),"]\n") } } ## pbeta(*, log=TRUE) {toms708} -- now improved tail behavior x <- c(.01, .10, .25, .40, .55, .71, .98) pbval <- c(-0.04605755624088, -0.3182809860569, -0.7503593555585, -1.241555830932, -1.851527837938, -2.76044482378, -8.149862739881) stopifnot(all.equal(pbeta(x, 0.8, 2, lower=FALSE, log=TRUE), pbval), all.equal(pbeta(1-x, 2, 0.8, log=TRUE), pbval)) qq <- 2^(0:1022) df.set <- c(0.1, 0.2, 0.5, 1, 1.2, 2.2, 5, 10, 20, 50, 100, 500) for(nu in df.set) { pqq <- pt(-qq, df = nu, log=TRUE) stopifnot(is.finite(pqq)) } ## PR#14230 -- more extreme beta cases {should no longer rely on denormalized} x <- (256:512)/1024 P <- pbeta(x, 3, 2200, lower.tail=FALSE, log.p=TRUE) stopifnot(is.finite(P), P < -600, -.001 < (D3P <- diff(P, diff = 3)), D3P < 0, diff(D3P) < 0) ## all but the first 43 where -Inf in R <= 2.9.1 stopifnot(All.eq(pt(2^-30, df=10), 0.50000000036238542)) ## = .5+ integrate(dt, 0,2^-30, df=10, rel.tol=1e-20) ## rbinom(*, size) gave NaN for large size up to R <= 2.6.1 M <- .Machine$integer.max set.seed(7) # as M is large, now "basically" rbinom(n, *) := qbinom(runif(n), *) : (tt <- table(rbinom(100, M, pr = 1e-9 )) ) # had values in {0,2} only (t2 <- table(rbinom(100, 10*M, pr = 1e-10)) ) stopifnot(0:6 %in% names(tt), sum(tt) == 100, sum(t2) == 100) ## no NaN there ## related qbinom() tests: (binomOk <- b64 && !(Lnx && usingMKL)) # not for MKL on RHEL {R-dev.: 2023-06-22} k <- 0:32 for(n in c((M+1)/2, M, 2*M, 10*M)) { for(pr in c(1e-8, 1e-9, 1e-10)) { nDup <- !duplicated( pb <- pbinom(k, n, pr) ) qb <- qbinom(pb[nDup], n, pr) pn1 <- pb[nDup] < if(binomOk) 1 else 1 - 3*.Machine$double.eps stopifnot(k[nDup][pn1] == qb[pn1]) ##^^^^^ fudge needed (Linux 32-b) } } ## qbinom() gave NaN in R 4.0.2 ## qf() with large df1, df2 and/or small p: x <- 0.01; f1 <- 1e60; f2 <- 1e90 stopifnot(qf(1/4, Inf, Inf) == 1, all.equal(1, 1e-18/ pf(qf(1e-18, 12,50), 12,50), tolerance = 1e-10), abs(x - qf(pf(x, f1,f2, log.p=TRUE), f1,f2, log.p=TRUE)) < 1e-4) ## qbeta(*, log.p) for "border" case: stopifnot(is.finite(q0 <- qbeta(-1e10, 50,40, log.p=TRUE)), 1 == qbeta(-1e10, 2, 3, log.p=TRUE, lower=FALSE)) ## infinite loop or NaN in R <= 2.7.0 ## phyper(x, 0,0,0), notably for huge x stopifnot(all(phyper(c(0:3, 1e67), 0,0,0) == 1)) ## practically infinite loop and NaN in R <= 2.7.1 (PR#11813) ## plnorm(<= 0, . , log.p=TRUE) stopifnot(plnorm(-1:0, lower.tail=FALSE, log.p=TRUE) == 0, plnorm(-1:0, lower.tail=TRUE, log.p=TRUE) == -Inf) ## was wrongly == 'log.p=FALSE' up to R <= 2.7.1 (PR#11867) ## pchisq(df=0) was wrong in 2.7.1; then, upto 2.10.1, P*(0,0) gave 1 stopifnot(pchisq(c(-1,0,1), df=0) == c(0,0,1), pchisq(c(-1,0,1), df=0, lower.tail=FALSE) == c(1,1,0), ## for ncp >= 80, gave values >= 1 in 2.10.0 pchisq(500:700, 1.01, ncp = 80) <= 1) ## dnbinom for extreme size and/or mu : mu <- 20 d <- dnbinom(17, mu=mu, size = 1e11*2^(1:10)) - dpois(17, lambda=mu) stopifnot(d < 0, diff(d) > 0, d[1] < 1e-10) ## was wrong up to 2.7.1 ## The fix to the above, for x = 0, had a new cancellation problem mu <- 1e12 * 2^(0:20) stopifnot(all.equal(1/(1+mu), dnbinom(0, size = 1, mu = mu), tolerance = 1e-13)) ## was wrong in 2.7.2 (only) mu <- sort(outer(1:7, 10^c(0:10,50*(1:6)))) NB <- dnbinom(5, size=1e305, mu=mu, log=TRUE) P <- dpois (5, mu, log=TRUE) stopifnot(abs(rErr(NB,P)) < 9*.Machine$double.eps)# seen 2.5* ## wrong in 3.1.0 and earlier ## Non-central F for large x x <- 1e16 * 1.1 ^ (0:20) dP <- diff(pf(x, df1=1, df2=1, ncp=20, lower.tail=FALSE, log=TRUE)) stopifnot(-0.047 < dP, dP < -0.0455) ## pf(*, log) jumped to -Inf prematurely in 2.8.0 and earlier ## Non-central Chi^2 density for large x stopifnot(0 == dchisq(c(Inf, 1e80, 1e50, 1e40), df=10, ncp=1)) ## did hang in 2.8.0 and earlier (PR#13309). ## qbinom() .. particularly for large sizes, small prob: p.s <- c(.01, .001, .1, .25) pr <- (2:20)*1e-7 sizes <- 1000*(5000 + c(0,6,16)) + 279 k.s <- 0:15; q.xct <- rep(k.s, each=length(pr)) for(sz in sizes) { for(p in p.s) { qb <- qbinom(p=p, size = sz, prob=pr) pb <- qpois (p=p, lambda = sz * pr) stopifnot(All.eq(qb, pb)) } pp.x <- outer(pr, k.s, function(pr, q) pbinom(q, size = sz, prob=pr)) qq.x <- apply(pp.x, 2, function(p) qbinom(p, size = sz, prob=pr)) stopifnot(qq.x == q.xct) } ## do_search() in qbinom() contained a thinko up to 2.9.0 (PR#13711) ## pbeta(x, a,b, log=TRUE) for small x and a is ~ log-linear x <- 2^-(200:10) for(a in c(1e-8, 1e-12, 16e-16, 4e-16)) for(b in c(0.6, 1, 2, 10)) { dp <- diff(pbeta(x, a, b, log=TRUE)) # constant approximately stopifnot(sd(dp) / mean(dp) < 0.0007) } ## had accidental cancellation '1 - w' ## qgamma(p, a) for small a and (hence) small p ## pgamma(x, a) for very very small a a <- 2^-seq(10,1000, .25) q.1c <- qgamma(1e-100,a,lower.tail=FALSE) q.3c <- qgamma(1e-300,a,lower.tail=FALSE) p.1c <- pgamma(q.1c[q.1c > 0], a[q.1c > 0], lower.tail=FALSE) p.3c <- pgamma(q.3c[q.3c > 0], a[q.3c > 0], lower.tail=FALSE) x <- 1+1e-7*c(-1,1); pg <- pgamma(x, shape = 2^-64, lower.tail=FALSE) stopifnot(qgamma(.99, .00001) == 0, abs(pg[2] - 1.18928249197237758088243e-20) < 1e-33, abs(diff(pg) + diff(x)*dgamma(1, 2^-64)) < 1e-13 * mean(pg), abs(1 - p.1c/1e-100) < 10e-13,# max = 2.243e-13 / 2.442 e-13 abs(1 - p.3c/1e-300) < 28e-13)# max = 7.057e-13 ## qgamma() was wrong here, orders of magnitude up to R 2.10.0 ## pgamma() had inaccuracies, e.g., ## pgamma(x, shape = 2^-64, lower.tail=FALSE) was discontinuous at x=1 stopifnot(all(qpois((0:8)/8, lambda=0) == 0)) ## gave Inf as p==1 was checked *before* lambda==0 ## extreme tail of non-central chisquare stopifnot(all.equal(pchisq(200, 4, ncp=.001, log.p=TRUE), -3.851e-42)) ## jumped to zero too early up to R 2.10.1 (PR#14216) ## left "extreme tail" lp <- pchisq(2^-(0:200), 100, 1, log=TRUE) stopifnot(is.finite(lp), lp < -184, all.equal(lp[201], -7115.10693158)) dlp <- diff(lp) dd <- abs(dlp[-(1:30)] - -34.65735902799) stopifnot(-34.66 < dlp, dlp < -34.41, dd < 1e-8)# 2.2e-10 64bit Lnx ## underflowed to -Inf much too early in R <= 3.1.0 for(e in c(0, 2e-16))# continuity at 80 (= branch point) stopifnot(all.equal(pchisq(1:2, 1.01, ncp = 80*(1-e), log=TRUE), c(-34.57369629, -31.31514671))) ## logit() == qlogit() on the right extreme: x <- c(10:80, 80 + 5*(1:24), 200 + 20*(1:25)) stopifnot(All.eq(x, qlogis(plogis(x, log.p=TRUE), log.p=TRUE))) ## qlogis() gave Inf much too early for R <= 2.12.1 ## Part 2: x <- c(x, seq(700, 800, by=10)) stopifnot(All.eq(x, qlogis(plogis(x, lower=FALSE, log.p=TRUE), lower=FALSE, log.p=TRUE))) # plogis() underflowed to -Inf too early for R <= 2.15.0 ## log upper tail pbeta(): x <- (25:50)/128 pbx <- pbeta(x, 1/2, 2200, lower.tail=FALSE, log.p=TRUE) d2p <- diff(dp <- diff(pbx)) b <- 2200*2^(0:50) y <- log(-pbeta(.28, 1/2, b, lower.tail=FALSE, log.p=TRUE)) stopifnot(-1094 < pbx, pbx < -481.66, -29 < dp, dp < -20, -.36 < d2p, d2p < -.2, all.equal(log(b), y+1.113, tolerance = .00002) ) ## pbx had two -Inf; y was all Inf for R <= 2.15.3; PR#15162 ## dnorm(x) for "large" |x| stopifnot(abs(1 - dnorm(35+3^-9)/ 3.933395747534971e-267) < 1e-15) ## has been losing up to 8 bit precision for R <= 3.0.x ## pbeta(x, ,, .., log): ldp <- diff(log(diff(pbeta(0.5, 2^-(90+ 1:25), 2^-60, log.p=TRUE)))) stopifnot(abs(ldp - log(1/2)) < 1e-9) ## pbeta(*, log) lost all precision here, for R <= 3.0.x (PR#15641) ## ## "stair function" effect (from denormalized numbers) a <- 43779; b <- 0.06728 x. <- .9833 + (0:100)*1e-6 px <- pbeta(x., a,b, log=TRUE) # plot(x., px) # -> "stair" d2. <- diff(dpx <- diff(px)) stopifnot(all.equal(px[1], -746.0986886924, tol=1e-12), 0.0445741 < dpx, dpx < 0.0445783, -4.2e-8 < d2., d2. < -4.18e-8) ## were way off in R <= 3.1.0 c0 <- system.time(p0 <- pbeta( .9999, 1e30, 1.001, log=TRUE)) cB <- max(.001, c0[[1]])# base time c1 <- system.time(p1 <- pbeta(1- 1e-9, 1e30, 1.001, log=TRUE)) c2 <- system.time(p2 <- pbeta(1-1e-12, 1e30, 1.001, log=TRUE)) stopifnot(all.equal(p0, -1.000050003333e26, tol=1e-10), all.equal(p1, -1e21, tol = 1e-6), all.equal(p2, -9.9997788e17), c(c1[[1]], c2[[1]]) < 1000*cB) ## (almost?) infinite loop in R <= 3.1.0 ## pbinom(), dbinom(), dhyper(),.. : R allows "almost integer" n for (FUN in c(function(n) dbinom(1,n,0.5), function(n) pbinom(1,n,0.5), function(n) dpois(n, n), function(n) dhyper(n+1, n+5,n+5, n))) try( lapply(sample(10000, size=1000), function(M) { ## invisible(lapply(sample(10000, size=1000), function(M) { n <- (M/100)*10^(2:20); if(anyNA(P <- FUN(n))) stop("NA for M=",M, "; 10ex=",paste((2:20)[is.na(P)], collapse=", "))})) ## check was too tight for large n in R <= 3.1.0 (PR#15734) ## [dpqr]beta(*, a,b) where a and/or b are Inf stopifnot(pbeta(.1, Inf, 40) == 0, pbeta(.5, 40, Inf) == 1, pbeta(.4, Inf,Inf) == 0, pbeta(.5, Inf,Inf) == 1, ## gave infinite loop (or NaN) in R <= 3.1.0 qbeta(.9, Inf, 100) == 1, # Inf.loop qbeta(.1, Inf, Inf) == 1/2)# NaN + Warning ## range check (in "close" cases): assertWarning(qN <- qbeta(2^-(10^(1:3)), 2,3, log.p=TRUE)) assertWarning(qn <- qbeta(c(-.1, -1e-300, 1.25), 2,3)) stopifnot(is.nan(qN), is.nan(qn)) ## lognormal boundary case sdlog = 0: p <- (0:8)/8; x <- 2^(-10:10) stopifnot(all.equal(qlnorm(p, meanlog=1:2, sdlog=0), qlnorm(p, meanlog=1:2, sdlog=1e-200)), dlnorm(x, sdlog=0) == ifelse(x == 1, Inf, 0)) ## qbeta(*, a,b) when a,b << 1 : can easily fail q1 <- qbeta(2^-28, 0.125, 2^-26) # gave 1000 Newton it + warning stopifnot(all.equal(2^-28, pbeta(q1, 0.125, 2^-26), tol= 2^-50)) a <- 1/8; b <- 2^-(4:200); alpha <- b/4 qq <- qbeta(alpha, a,b)# gave warnings intermediately pp <- pbeta(qq, a,b) stopifnot(pp > 0, diff(pp) < 0, ## pbeta(qbeta(alpha,*),*) == alpha: abs(1 - pp/alpha) < 4e-15)# seeing 2.2e-16 ## orig. qbeta() using *many* Newton steps; case where we "know the truth" a <- 25; b <- 6; x <- 2^-c(3:15, 100, 200, 250, 300+100*(0:7)) pb <- c(## via Rmpfr's roundMpfr(pbetaI(x, a,b, log.p=TRUE, precBits = 2048), 64) : -40.7588797271766572448, -57.7574063441183625303, -74.9287878018119846216, -92.1806244636893542185, -109.471318248524419364, -126.781111923947395655, -144.100375042814531426, -161.424352961544612370, -178.750683324909148353, -196.078188674895169383, -213.406281209657976525, -230.734667259724367416, -248.063200048177428608, -1721.00081201679567511, -3453.86876341665894863, -4320.30273911659058550, -5186.73671481652222237, -6919.60466621638549567, -8652.47261761624876897, -10385.3405690161120427, -12118.2085204159753165, -13851.0764718158385902, -15583.9444232157018631, -17316.8123746155651368) stopifnot(all.equal(pb, pbeta(x,a,b, log.p=TRUE), tol=8e-16))# seeing {1.5|1.6|2.0}e-16 qp <- qbeta(pb, a,b, log.p=TRUE) ## x == qbeta(pbeta(x, *), *) : stopifnot(qp > 0, all.equal(x, qp, tol= 1e-15))# seeing {2.4|3.3}e-16 ## qbeta(), PR#15755 a1 <- 0.0672788; b1 <- 226390 curve(pbeta(x, a1,b1), from=1e-50, n=4001, log="x", main = sprintf("(a1=%g, b1=%g) -- needs log-x scale",a1,b1)) ## zoom into part where it jumps from const = 0 to const = 5.562685e-309 to regular growth: curve(qbeta(x, a1,b1), from=1e-21, to=6e-21, type="o", cex=1/4) -> qb1 rle(qb1$y) p <- 0.6948886 qp <- qbeta(p, a1,b1) stopifnot(qp < 2e-8, # was '1' (with a warning) in R <= 3.1.0 All.eq(p, pbeta(qp, a1,b1))) ## less extreme example, same phenomenon: a <- 43779; b <- 0.06728 stopifnot(All.eq(0.695, pbeta(qbeta(0.695, b,a), b,a))) x <- -exp(seq(0, 14, by=2^-9)) qx <- qbeta(x, a,b, log.p=TRUE)# used to be slow pqx <- pbeta(qx, a,b, log=TRUE) stopifnot(diff(qx) < 0, all.equal(x, pqx, tol= 2e-15)) # seeing {3.51|3.54}e-16 ## note that qx[x > -exp(2)] is too close to 1 to get full accuracy: i2 <- x > -exp(2); all.equal(x[i2], pqx[i2], tol= 0)#-> 5.849e-12 ## was Inf, and much slower, for R <= 3.1.0 x3 <- -(15450:15700)/2^11 pq3 <- pbeta(qbeta(x3, a,b, log.p=TRUE), a,b, log=TRUE) stopifnot(mean(abs(pq3-x3)) < 4e-12,# 1.46e-12 max (abs(pq3-x3)) < 8e-12)# 2.95e-12 ## .a <- .2; .b <- .03; lp <- -(10^-(1:323)) qq <- qbeta(lp, .a,.b, log=TRUE) # warnings in R <= 3.1.0 assertWarning(qN <- qbeta(.5, 2,3, log.p=TRUE)) assertWarning(qn <- qbeta(c(-.1, 1.25), 2,3)) stopifnot(1-qq < 1e-15, is.nan(qN), is.nan(qn))# typically qq == 1 exactly ## failed in intermediate versions ## a <- 2^-8; b <- 2^(200:500) pq <- pbeta(qbeta(1/8, a, b), a, b) stopifnot(abs(pq - 1/8) < 1/8) ## whereas qbeta() would underflow to 0 "too early", for R <= 3.1.0 # ## very extreme tails on both sides x <- c(1e-300, 1e-12, 1e-5, 0.1, 0.21, 0.3) stopifnot(0 == qbeta(x, 2^-12, 2^-10))## gave warnings a <- 10^-(8:323) length(uq <- unique(qb <- qbeta(0.95, a, 20)))# typically just 1 ! head(uq) # 5.562685e-309 ## had warnings and wrong value +1; also NaN ct2 <- system.time(q2 <- qbeta(0.95, a,a))[1] stopifnot(is.finite(qb), qb < 1e-300, q2 == 1) if(ct2 > 0.020) { cat("system.time:\n"); print(ct2) } ## had warnings and was much slower for R <= 3.1.0 ## qt(p, df= Inf, ncp) <==> qnorm(p, m=ncp) p <- (0:32)/32 stopifnot(all.equal(qt(p, df=Inf, ncp=5), qnorm(p, m=5))) ## qt(*, df=Inf, .) gave NaN in R <= 3.2.1 ## rhyper(*, ); PR#16489 ct3 <- system.time(N <- rhyper(100, 8000, 1e9-8000, 1e6))[1] table(N) stopifnot(exprs = { identical(c(table(N)), `names<-`(as.integer(c(1, 1, 2, 2, 8, 21, 9, 18, 15, 10, 1, 7, 1, 2, 2)), as.character((0:15)[-2]))) abs(mean(N) - 8) < 1.5 }) if(ct3 > 0.02) { cat("system.time:\n"); print(ct3) } ## N were all 0 and took very long for R <= 3.2.1 set.seed(17) stopifnot(rhyper(1, 3024, 27466, 251) == 25, rhyper(1, 329, 3059, 225) == 22) ## failed for a day after a "thinko" in the above bug fix. ## *chisq(*, df=0, ncp=0) == Point mass at 0 stopifnot(rchisq(32, df=0, ncp=0) == 0, dchisq((0:16)/16, df=0, ncp=0) == c(Inf, rep(0, 16))) ## gave all NaN's for R <= 3.2.2 ## pchisq(*, df=0, ncp > 0, log.p=TRUE) : th <- 10*c(1:9,10^c(1:3,7)) pp <- pchisq(0, df = 0, ncp=th, log.p=TRUE) stopifnot(all.equal(pp, -th/2, tol=1e-15)) ## underflowed at about th ~= 60 in R <= 3.2.2 ## pnbinom (-> C's bratio()) op <- options(warn = 1)# -- NaN's giving warnings L <- 1e308; p <- suppressWarnings(pnbinom(L, L, mu = 5)) # NaN or 1 (for 64 / 32 bit) stopifnot(is.nan(p) || p == 1) ## gave infinite loop on some 64b platforms in R <= 3.2.3 ## [dpqr]nbinom(*, mu, size=Inf) -- PR#16727 ## (extended wrt changes for PR#15628 ) L <- 1e308; mu <- 5 x <- c(0:3, 1e10, 1e100, L, Inf) MxM <- .Machine$double.xmax xL <- c(2^(1000+ c(1:23, 23.5, 23.9, 24-1e-13)), MxM, Inf) (dP <- dpois(xL, mu)) # all 0 (pP <- ppois(xL, mu)) # all 1 stopifnot(dP == 0, pP == 1, identical(pP, pgamma(mu, xL + 1, 1., lower.tail=FALSE))) ## cbind(xL, dP, pP) (dLmI <- dnbinom(xL, mu = 1, size = Inf)) # all == 0 ## FIXME ?!: MxM/2 seems +- ok ?? (dLmM <- dnbinom(xL, mu = 1, size = MxM)) # all NaN but the last (dLpI <- dnbinom(xL, prob=1/2, size = Inf))# ditto (dLpM <- dnbinom(xL, prob=1/2, size = MxM))# ditto d <- dnbinom(x, mu = mu, size = Inf) # gave NaN (for 0 and L), now all 0 p <- pnbinom(x, mu = mu, size = Inf) # gave all NaN, now uses ppois(x, mu) pp <- (0:16)/16 q <- qnbinom(pp, mu = mu, size = Inf) # gave all NaN set.seed(1); NI <- rnbinom(32, mu = mu, size = Inf)# gave all NaN set.seed(1); N2 <- rnbinom(32, mu = mu, size = L ) stopifnot(exprs = { all.equal(d, c(0.006737947, 0.033689735, 0.0842243375, 0.140373896, 0,0,0,0), tol = 9e-9)# 7.6e-10 all.equal(p, c(0.006737947, 0.040427682, 0.1246520195, 0.265025915, 1,1,1,1), tol = 9e-9)# 7.3e-10 all.equal(d, dpois(x, mu))# current implementation: even identical() all.equal(p, ppois(x, mu)) q == c(0, 2, 3, 3, 3, 4, 4, 4, 5, 5, 6, 6, 6, 7, 8, 9, Inf) q == qpois(pp, mu) identical(NI, N2) }) options(op) ## size = Inf -- mostly gave NaN in R <= 3.2.3 ## qpois(p, *) for invalid 'p' should give NaN -- PR#16972 stopifnot(is.nan(suppressWarnings(c(qpois(c(-2,3, NaN), 3), qpois(1, 3, log.p=TRUE), qpois(.5, 0, log.p=TRUE), qpois(c(-1,pi), 0))))) ## those in the 2nd line gave 0 in R <= 3.3.1 ## Similar but different for qgeom(): stopifnot(qgeom((0:8)/8, prob=1) == 0, ## p=1 gave Inf in R <= 3.3.1 is.nan(suppressWarnings(qgeom(c(-1/4, 1.1), prob=1)))) ## all our RNG r() functions: ##' catch all: value and warnings or error <-- demo(error.catching) : tryCatch.W.E <- function(expr) { W <- NULL w.handler <- function(w){ # warning handler W <<- w invokeRestart("muffleWarning") } list(value = withCallingHandlers(tryCatch(expr, error = function(e) e), warning = w.handler), warning = W) } .stat.ns <- asNamespace("stats") Ns <- 4 for(dist in PDQR) { fn <- paste0("r",dist) cat(sprintf("%-9s(%d, ..): ", fn, Ns)) F <- get(fn, envir = .stat.ns) nArg <- length(fms <- formals(F)) if(dist %in% c("nbinom", "gamma")) ## cannot specify *both* 'prob' & 'mu' / 'rate' & 'scale' nArg <- nArg - 1 nA1 <- nArg - 1 # those beside the first (= 'n' mostly) expected <- rep(if(dist %in% PDQRinteg) NA_integer_ else NaN, Ns) for(ia in seq_len(nA1)) { aa <- rep(list(1), nA1) aa[[ia]] <- NA cat(ia,"") R <- tryCatch.W.E( do.call(F, c(Ns, aa)) ) if(!inherits(R$warning, "simpleWarning")) stop(" .. did *NOT* give a warning! ") if(!(identical(R$value, expected))) { ## allow NA/NaN mismatch in these cases for now: if(!(dist %in% c("beta","f","t") && all(is.na(R$value)))) stop(" .. not giving expected NA/NaN's ") } } cat(" [Ok]\n") } ## qbeta() in very asymmetric cases sh2 <- 2^seq(9,16, by=1/16) p <- 1e-10 qbet <- qbeta(p, 1.5, shape2=sh2, lower.tail=FALSE) plot(sh2, pbeta(qbet, 1.5, sh2, lower.tail=FALSE)/p -1 -> rE, log="x", main="rel.Error") dqb <- diff(qbet); d2qb <- diff(dqb); d3qb <- diff(d2qb) stopifnot(all.equal(qbet[[1]], 0.047206901483498, tol=1e-12), print(max(abs(rE))) < 1e-12, # Lx 64b: 2.4e-13 0 > dqb, dqb > -0.002, 0 < d2qb, d2qb < 0.00427, -3.2e-8 > d3qb, d3qb > -3.1e-6, diff(d3qb) > 1e-9) ## had discontinuity (from wrong jump out of Newton) in R <= 3.3.2 ## rt() [PR#17306]; rf() and rbeta() [PR#17375] with non-scalar 'ncp' nc <- c(NA, 1); iN <- is.na(rep_len(nc, 3)) ## each gives warning "NAs produced": assertWarning(T <- rt (3, 4, ncp = nc)) assertWarning(F <- rf (3, 4,5, ncp = nc)) assertWarning(B <- rbeta(3, 4,5, ncp = nc)) stopifnot(identical(iN, is.na(T)), identical(iN, is.na(F)), identical(iN, is.na(B))) ## was not handled correctly, notably with NA's in ncp, in R <= 3.4.(2|3) ## check old version of walker_Probsample is being used for old sample kind suppressWarnings(RNGversion("3.5.0")) set.seed(12345) p <- c(2, rep(1, 200)) x <- sample(length(p), 100000, prob = p, replace = TRUE) stopifnot(sum(x == 1) == 994) ## check for failure of new walker_Probsample RNGversion("3.6.0") set.seed(12345) epsilon <- 1e-10 p201 <- proportions( rep( c(1, epsilon), c(201, 999-201))) x <- sample(length(p201), 100000, prob = p201, replace = TRUE) stopifnot(sum(x <= 201) == 100000) ## had if(!(onWindows && arch == "x86")) ## PR#17577 - dgamma(x, shape) for shape < 1 (=> +Inf at x=0) and very small x stopifnot(exprs = { all.equal(dgamma(2^-1027, shape = .99 , log=TRUE), 7.1127667376, tol=1e-10) all.equal(dgamma(2^-1031, shape = 1e-2, log=TRUE), 702.8889158, tol=1e-10) all.equal(dgamma(2^-1048, shape = 1e-7, log=TRUE), 710.30007699, tol=1e-10) all.equal(dgamma(2^-1048, shape = 1e-7, scale = 1e-315, log=TRUE), 709.96858768, tol=1e-10) }) ## all gave Inf in R <= 3.6.1 ## } else cat("PR#17577 bug fix not checked, as it may not work on this platform\n") ## if(!(onWindows && arch == "x86")) { ## This gave a practically infinite loop (on 64-bit Lnx, Windows; not in 32-bit) suppressWarnings(# typically warns, but not on Apple clang 14.0.3 p <- pchisq(1.00000012e200, df=1e200, ncp=100) ) stopifnot(p == 1) ## } ## Extreme tails for qnorm(*, log.p=TRUE) : qs <- 2^seq(0, 155, by=1/8) lp <- pnorm( qs, log.p=TRUE, lower.tail=FALSE) lp. <- pnorm(-qs, log.p=TRUE) stopifnot(all.equal(lp, lp., tol= 4e-16)) # actually exactly via code-identity ## Both these gave NaNs (and warned about it): qpU <- qnorm(lp, log.p=TRUE, lower.tail=FALSE) qp. <- qnorm(lp, log.p=TRUE) ## Show the (mostly) small differences : all.equal( qs, qpU, tol=0) all.equal(-qs, qp., tol=0) all.equal(-qp.,qpU, tol=0) # typically TRUE (<==> exact equality) ## however, range(qpU/qs - 1) # -5.68e-6 5.41e-6 in R <= 4.2.1 stopifnot(exprs = { all.equal( qs, qpU, tol=1e-15) all.equal(-qs, qp., tol=1e-15) all.equal(-qp., qpU, tol=1e-15)# diff of 4.71e-16 in 4.1.0 w/icc (Eric Weese) max(abs(qpU/qs - 1)) < 1e-15 # see 4.44e-16 {was 5.68e-6 in R <= 4.2.1; much larger in R <= 4.0.x) }) ## both failed very badly in R <= 4.0.x ## pnorm(x, log.p=TRUE) now works for larger x, |x| <= 1.896..e154 {before: |x| < 1.34e154} x <- seq(134.5, 189, by=.5) px <- pnorm(-x * 1e152, log.p=TRUE)# all these underflowed previously stopifnot(exprs = { all.equal(-1.79769313486073e+308, pnorm(-1.896150381621e154, log.p=TRUE), tol=1e-14) is.finite(px) abs(1 - diff(diff(px)) / -2.5e303) < 3e-11 * (1 + (.Machine$sizeof.longdouble < 12)) }) ## all these where -Inf in R <= 4.0.x ## qnbinom(*, size=, mu=) -- PR#18095: qi <- 0:16 qiL <- c(0:99, 10L*(10:20), 50L*(5:20)) verbose <- interactive() for(N in c(10^c(-300, -100, -20, -9:-7, -2, 2, 7:9), 9.1e15, 10^c(20,30,50, 100, 300))) { if(verbose) cat("N =", N,": ") pbi <- pnbinom(qi, size=N, mu=1) pbiU <- pnbinom(qiL, size=N, mu=1, lower.tail=FALSE, log.p=TRUE) stopifnot(exprs = { ## whne pbinom(.) gave 1, quntile is Inf pbi == 1 | qi == qnbinom(pbi, size=N, mu=1) qiL == qnbinom(pbiU, size=N, mu=1, lower.tail=FALSE, log.p=TRUE) }) if(verbose) cat("[Ok]\n") } ## qnbinom(*) gave all 0 in R <= 4.1.0 ## ## ## Fix qnbinom(*, scale = <>), partly reported to R-devel, Aug.7, 2020 ## by Constantin Ahlmann-Eltze: size=1e-10 already took 30 sec on his computer: st <- system.time(qn <- qnbinom(0.5, mu = 3, size = 2^-(10:59)))[[1]] stopifnot(exprs = { st < .5 # = 500 ms; observed 0 sec, i.e., '< 1 ms' qn == 0 print(qnbinom(.9999, mu=3, size=1e-4)) == 7942 # was 8044 (MM on R-devel) print(qnbinom(.9999, mu=1e5, size=1e10)) == 101178 # was off by 1 { ## Even larger: ^^^^ large size size <- 1e13; mu <- size/4; k <- 2500002265485 + (-100:100) pnb <- pnbinom( k , mu=mu, size=size) # around 0.9 qnb <- qnbinom(pnb, mu=mu, size=size) qnb == k ## in R <= 4.0.2, qnb == k+23 wrongly } { ## small size, large k : needs "p fuzz factor 16 instead of 64 : k <- 75000 + -1000:1000; size <- 1e-9; mu <- 30 pnb <- pnbinom( k , mu=mu, size=size) # = 0.999999987678 (for k=75'000) qnb <- qnbinom(pnb, mu=mu, size=size) qnb == k ## in R <= 4.0.2, qnb == k-1 wrongly } }) ## For log.p=TRUE and or lower.tail=FALSE ## ## all three qpois(), qbinom() & qnbinom() were *far* from good in R <= 4.1.0 ## PR#18072 : dnbinom() underflow stopifnot(dnbinom(1:40, size=2^58, prob = 1) == 0) ## gave mostly 1 in R <= 4.1.0 x <- unique(sort(c(1:12, 15, outer(c(1,2,5), 10^(1:11))))) sz <- 2^70 ; prb <- .9999999 summary(dn <- dnbinom(x, size=sz, prob = prb, log=TRUE)) dL <- 118059167912526.5 summary(dl.dn1 <- diff(log(dn[-1] + dL))) stopifnot(dn + dL > 0, 0.09 < dl.dn1, dl.dn1 < 0.93) ## accuracy loss of 6 and more digits in R <= 4.1.0 ##---- reverse case, very *small* size --------------- dS <- dnbinom(1:90, size=1e-15, mu=200, log=TRUE) d4S <- diff(d3S <- diff(ddS <- diff(dS))) stopifnot(-39.1 < dS, dS < -34.53 , -0.7 < ddS, ddS < -0.01116 , 0.000126 < d3S, d3S < 0.287683 , -0.17 < d4S, d4S < -2.8859e-6 , all.equal(c(-48.40172, -49.155492, -49.905797, -50.653012, -51.397452), dnbinom(16:20, size=1e-15, prob=1/2, log=TRUE)) ) ## failed in R 4.1.1 (and R-devel) only ## PR#15628 ebd0() for dpois() [for now] \\ dpois(x, *) for 2\pi x >= DMAX stopifnot(exprs = { print(abs(dpois( 40, 7.5)/ 6.81705586862060455e-17 -1)) < 6e-16 # 2.66e-15 in R <= 4.1.1 print(abs(dpois(1400,1000)/ 1.46677334419656833e-33 -1)) < 8e-16 # 1.71e-14 in R <= 4.1.1 }) ## very large x: x <- trunc(2^(1000+ head(seq(1,24, by=1/64), -1))) L <- tail(x,1) dpxx <- dpois(x,x) ## had ended in many 0's ldpxx <- dpois(x,x, log=TRUE) # ... -Inf (d <- mean(dlp <- diff(ldpxx)))# -0.005415 stopifnot(exprs = { dpxx > 0 is.finite(ldpxx) print(abs(print(dpois(L,L))/ (1/sqrt(2*pi)/sqrt(L)) -1)) < 1e-15 # see 1.11e-16 abs(range(dlp) - d) < 1e-12 # seen 4.4e-14, was NaN in R <= 4.1.1 all.equal(ldpxx, log(dpxx), tol = 1e-15) }) ## dpois(x,x) underflowed to zero in R <= 4.1.1 for such large x. ## PR#18642 -- dgeom() accuracy --> improved via dbinom_raw(x, n, prob) for x=0 and x=n x <- c(159, 171, 183, 201) tru1 <- c(4.64906012307645596e-240, 4.03241686836417614e-258, 3.4975641032381073e-276, 2.82530978257810403e-303) (xs <- 44400 + sort(c(outer(c(10,26), 29*(0:2), `+`)))) tru2 <- c(2.850864888117265e-306, 2.2158779845990397e-306, 1.8056489670203573e-306, 1.4034680530148749e-306, 1.1436417789181365e-306, 8.8891292278883415e-307) stopifnot(exprs = { print(abs(1 - dgeom(x, 31/32) / tru1) * 2^52) <= 4 # see 0 0 0 0 (Linux F 38, gcc) print(abs(1 - dgeom(xs, 1/64) / tru2) * 2^52) <= if(onWindows) 800 else 4 ## on Windows: 406.5 408.0 407.0 407.0 408.0 407.5; see 0 0 0 0 0 0 ( " ) ## Reprex for Windows: print(2^52 * abs(1 - dgeom(44410, 1/64) / 2.850864888117265e-306)) < 606 # -> 406.5; R 4.3.2 Lnx: 252 }) # in R <= 4.3.z, relErr * 2^52 were (238 246 254 10) and (252 242.5 252 242 252 242) ## PR#18672 -- x = 0|1 when one or both shape = 0 stopifnot(exprs = { pbeta(0, 0, 3) == 1 # gave 0 pbeta(1, .1, 0) == 1 # gave 0 pbeta(1.1, 3,0) == 1 # gave 0 pbeta(0, 0, 0) == 0.5 # gave 0, should give 0.5 pbeta(1, 0, 0) == 1 # gave 0.5 }) ## PR#18640 -- stirlerr(x) concerns (for *non* half-integer x) -- visible in dgamma() sh <- 465/32 # = 14.53125 x0 <- 1/4 + 8:20 dg1 <- dgamma(x0, sh) ## 'TRUE' values { = dput(asNumeric(Rmpfr::dgamma(mpfr(x0, 512), sh)), control="digits17") } : dgM <- c(0.026214161736344995, 0.045350212095058476, 0.066917055544970391, 0.086754619023375584, 0.10102573200716865, 0.10746812620489894, 0.10581786016571779, 0.097460173977057932, 0.084678318901885929, 0.069890902339799707, 0.055117163281675971, 0.041733282935808788, 0.030464801624510086) relE <- dg1/dgM - 1 relE * 2^53 # 2 2 -2 0 0 2 -1 0 0 2 0 0 2 // was in {-95 : -91} in R <= 4.3.* stopifnot(abs(relE) < 9e-16) # max{x86_64}: 2.22e-16 ## PR#18711 -- qbinom() - inversion of pbinom() ## but probably also fixing qpois(), qnbinom() cases sz <- 6040:6045 prb <- 0.995 (qb6 <- qbinom(p = 0.05, size = sz, prob = prb)) (pqb6 <- pbinom(qb6, size = sz, prob = prb)) (pqb6_1 <- pbinom(qb6-1, size = sz, prob = prb)) stopifnot(exprs = { qb6 == c(6001:6004,6004:6005) # not in R 4.4.0, nor 4.1.1 1 > pqb6 & pqb6 >= 0.05 # " 0.05 > pqb6_1 & pqb6_1 >= 0.035# " }) ## was wrong for R versions in [4.1.1, 4.4.0] cat("Time elapsed: ", proc.time() - .ptime,"\n")