scipy.special.nbdtr#
- scipy.special.nbdtr(k, n, p, out=None) = <ufunc 'nbdtr'>#
Negative binomial cumulative distribution function.
Returns the sum of the terms 0 through k of the negative binomial distribution probability mass function,
\[F = \sum_{j=0}^k {{n + j - 1}\choose{j}} p^n (1 - p)^j.\]In a sequence of Bernoulli trials with individual success probabilities p, this is the probability that k or fewer failures precede the nth success.
- Parameters:
- karray_like
The maximum number of allowed failures (nonnegative int).
- narray_like
The target number of successes (positive int).
- parray_like
Probability of success in a single event (float).
- outndarray, optional
Optional output array for the function results
- Returns:
- Fscalar or ndarray
The probability of k or fewer failures before n successes in a sequence of events with individual success probability p.
See also
nbdtrc
Negative binomial survival function
nbdtrik
Negative binomial quantile function
scipy.stats.nbinom
Negative binomial distribution
Notes
If floating point values are passed for k or n, they will be truncated to integers.
The terms are not summed directly; instead the regularized incomplete beta function is employed, according to the formula,
\[\mathrm{nbdtr}(k, n, p) = I_{p}(n, k + 1).\]Wrapper for the Cephes [1] routine
nbdtr
.The negative binomial distribution is also available as
scipy.stats.nbinom
. Usingnbdtr
directly can improve performance compared to thecdf
method ofscipy.stats.nbinom
(see last example).References
[1]Cephes Mathematical Functions Library, http://www.netlib.org/cephes/
Examples
Compute the function for
k=10
andn=5
atp=0.5
.>>> import numpy as np >>> from scipy.special import nbdtr >>> nbdtr(10, 5, 0.5) 0.940765380859375
Compute the function for
n=10
andp=0.5
at several points by providing a NumPy array or list for k.>>> nbdtr([5, 10, 15], 10, 0.5) array([0.15087891, 0.58809853, 0.88523853])
Plot the function for four different parameter sets.
>>> import matplotlib.pyplot as plt >>> k = np.arange(130) >>> n_parameters = [20, 20, 20, 80] >>> p_parameters = [0.2, 0.5, 0.8, 0.5] >>> linestyles = ['solid', 'dashed', 'dotted', 'dashdot'] >>> parameters_list = list(zip(p_parameters, n_parameters, ... linestyles)) >>> fig, ax = plt.subplots(figsize=(8, 8)) >>> for parameter_set in parameters_list: ... p, n, style = parameter_set ... nbdtr_vals = nbdtr(k, n, p) ... ax.plot(k, nbdtr_vals, label=rf"$n={n},\, p={p}$", ... ls=style) >>> ax.legend() >>> ax.set_xlabel("$k$") >>> ax.set_title("Negative binomial cumulative distribution function") >>> plt.show()
The negative binomial distribution is also available as
scipy.stats.nbinom
. Usingnbdtr
directly can be much faster than calling thecdf
method ofscipy.stats.nbinom
, especially for small arrays or individual values. To get the same results one must use the following parametrization:nbinom(n, p).cdf(k)=nbdtr(k, n, p)
.>>> from scipy.stats import nbinom >>> k, n, p = 5, 3, 0.5 >>> nbdtr_res = nbdtr(k, n, p) # this will often be faster than below >>> stats_res = nbinom(n, p).cdf(k) >>> stats_res, nbdtr_res # test that results are equal (0.85546875, 0.85546875)
nbdtr
can evaluate different parameter sets by providing arrays with shapes compatible for broadcasting for k, n and p. Here we compute the function for three different k at four locations p, resulting in a 3x4 array.>>> k = np.array([[5], [10], [15]]) >>> p = np.array([0.3, 0.5, 0.7, 0.9]) >>> k.shape, p.shape ((3, 1), (4,))
>>> nbdtr(k, 5, p) array([[0.15026833, 0.62304687, 0.95265101, 0.9998531 ], [0.48450894, 0.94076538, 0.99932777, 0.99999999], [0.76249222, 0.99409103, 0.99999445, 1. ]])