Statistics: random number generators, PDF and CDF functions, and hypothesis tests. The module also includes some basic statistical functions such as mean, variance, skew, and etc.
concordant x y calculates the number of concordant pairs in the given arrays x and y. A pair of observations \((x_i, y_i)\) and \((x_j, y_j)\) is concordant if both elements in one pair are either greater than or less than both elements in the other pair, i.e., \((x_i > x_j) \land (y_i > y_j)\) or \((x_i < x_j) \land (y_i < y_j)\).
discordant x y calculates the number of discordant pairs in the given arrays x and y. A pair of observations \((x_i, y_i)\) and \((x_j, y_j)\) is discordant if one element in one pair is greater than the corresponding element in the other pair, and the other element is smaller, i.e., \((x_i > x_j) \land (y_i < y_j)\) or \((x_i < x_j) \land (y_i > y_j)\).
quantile x p returns the p quantile of the data x. x does not need to be sorted beforehand. When computing several quantiles on the same data, it is more efficient to "pre-apply" the function, as in ``let f = quantile x in List.map f 0.1 ; 0.5 ``, in which case the data is sorted only once.
Raises Invalid_argument if p is not between 0 and 1.
sort x sorts the elements in the x in increasing order if inc = true, otherwise in decreasing order if inc=false. By default, inc is true. Note a copy is returned, the original data is not modified.
The ranking order is from the smallest one to the largest. For example rank [|54.; 74.; 55.; 86.; 56.|] returns [|1.; 4.; 2.; 5.; 3.|]. Note that the ranking starts with one!
ties_strategy controls which ranks are assigned to equal values:
Average the mean of ranks should be assigned to each value. Default.
Min the minimum of ranks is assigned to each value.
Max the maximum of ranks is assigned to each value.
Type for computed histograms, with optional weighted counts and normalized counts.
Sourceval histogram :
[ `Bins of float array| `N of int ]->?weights:float array->float array->histogram
histogram bins x creates a histogram from values in x. If bins matches `N n it will construct n equally spaced bins from the minimum to the maximum in x. If bins matches `Bins b, b is taken as the sorted array of boundaries of adjacent bin intervals. Bin boundaries are taken as left-inclusive, right-exclusive, except for the last bin which is also right-inclusive. Values outside the bins are dropped silently.
histogram bins ~weights x creates a weighted histogram with the given weights which must match x in length. The bare counts are also provided.
Returns a histogram including the n+1 bin boundaries, n counts and weighted counts if applicable, but without normalisation.
Sourceval histogram_sorted :
[ `Bins of float array| `N of int ]->?weights:float array->float array->histogram
histogram_sorted bins x is like histogram but assumes that x is sorted already. This increases efficiency if there are less bins than data. Undefined results if x is not in fact sorted.
normalize hist calculates a probability mass function using hist.weighted_counts if present, otherwise using hist.counts. The result is stored in the normalised_counts field and sums to one.
normalize_density hist calculates a probability density function using hist.weighted_counts if present, otherwise using hist.counts. The result is normalized as a density that is piecewise constant over the bin intervals. That is, the sum over density times corresponding bin width is one. If bins are infinitely wide, their density is 0 and the sum over width times density of all finite bins is the total weight in the finite bins. The result is stored in the density field.
ecdf x returns (x',f) which are the empirical cumulative distribution function f of x at points x'. x' is just x sorted in increasing order with duplicates removed. The function does not support nan values in the array x.
normalise_pdf arr takes an array of floats arr representing a probability density function (PDF) and returns a new array where the values are normalized so that the sum of the array equals 1.
parameterarr
The input array representing the unnormalized PDF.
returns
A new array of the same length as arr, with values normalized to ensure the sum is 1.
tukey_fences ?k x returns a tuple of the lower and upper boundaries for values that are not outliers. k defaults to the standard coefficient of 1.5. For first and third quartiles Q1 and `Q3`, the range is computed as follows:
gaussian_kde x is a Gaussian kernel density estimator. The estimation of the pdf runs in `O(sample_size * n_points)`, and returns an array tuple (a, b) where a is a uniformly spaced points from the sample range at which the density function was estimated, and b is the estimates at these points.
Bandwidth selection rules is as follows: * Silverman: use `rule-of-thumb` for choosing the bandwidth. It defaults to 0.9 * min(SD, IQR / 1.34) * n^-0.2. * Scott: same as Silverman, but with a factor, equal to 1.06.
metropolis_hastings target_density initial_state num_samples performs the Metropolis-Hastings algorithm to generate samples from a target distribution.
The Metropolis-Hastings algorithm is a Markov Chain Monte Carlo (MCMC) method that generates a sequence of samples from a probability distribution for which direct sampling is difficult. The algorithm uses a proposal distribution to explore the state space, accepting or rejecting proposed moves based on the ratio of the target densities.
parametertarget_density
A function that computes the density (up to a normalizing constant) of the target distribution at a given point.
parameterinitial_state
The starting point of the Markov chain, represented as an array of floats.
parameternum_samples
The number of samples to generate.
returns
A 2D array where each row represents a sample generated by the Metropolis-Hastings algorithm.
gibbs_sampling conditional_sampler initial_state num_samples performs Gibbs sampling to generate samples from a multivariate distribution.
Gibbs sampling is a Markov Chain Monte Carlo (MCMC) algorithm used to generate samples from a joint distribution when direct sampling is difficult. It works by iteratively sampling each variable from its conditional distribution, given the current values of all other variables.
parameterconditional_sampler
A function that takes the current state of the variables (as a float array) and the index of the variable to sample, and returns a new value for that variable sampled from its conditional distribution.
parameterinitial_state
The starting point of the Markov chain, represented as an array of floats.
parameternum_samples
The number of samples to generate.
returns
A 2D array where each row represents a sample generated by the Gibbs sampling algorithm.
z_test ~mu ~sigma ~alpha ~side x returns a test decision for the null hypothesis that the data x comes from a normal distribution with mean mu and a standard deviation sigma, using the z-test of alpha significance level. The alternative hypothesis is that the mean is not mu.
The result (h,p,z) : h is true if the test rejects the null hypothesis at the alpha significance level, and false otherwise. p is the p-value and z is the z-score.
t_test ~mu ~alpha ~side x returns a test decision of one-sample t-test which is a parametric test of the location parameter when the population standard deviation is unknown. mu is population mean, alpha is the significance level.
t_test_paired ~alpha ~side x y returns a test decision for the null hypothesis that the data in x – y comes from a normal distribution with mean equal to zero and unknown variance, using the paired-sample t-test.
t_test_unpaired ~alpha ~side ~equal_var x y returns a test decision for the null hypothesis that the data in vectors x and y comes from independent random samples from normal distributions with equal means and equal but unknown variances, using the two-sample t-test. The alternative hypothesis is that the data in x and y comes from populations with unequal means.
equal_var indicates whether two samples have the same variance. If the two variances are not the same, the test is referred to as Welche's t-test.
ks_test ~alpha x f returns a test decision for the null hypothesis that the data in vector x comes from independent random samples of the distribution with CDF f. The alternative hypothesis is that the data in x comes from a different distribution.
The result (h,p,d) : h is true if the test rejects the null hypothesis at the alpha significance level, and false otherwise. p is the p-value and d is the Kolmogorov-Smirnov test statistic.
ks2_test ~alpha x y returns a test decision for the null hypothesis that the data in vectors x and y come from independent random samples of the same distribution. The alternative hypothesis is that the data in x and y are sampled from different distributions.
The result (h,p,d): h is true if the test rejects the null hypothesis at the alpha significance level, and false otherwise. p is the p-value and d is the Kolmogorov-Smirnov test statistic.
var_test ~alpha ~side ~variance x returns a test decision for the null hypothesis that the data in x comes from a normal distribution with input variance, using the chi-square variance test. The alternative hypothesis is that x comes from a normal distribution with a different variance.
jb_test ~alpha x returns a test decision for the null hypothesis that the data x comes from a normal distribution with an unknown mean and variance, using the Jarque-Bera test.
fisher_test ~alpha ~side a b c d fisher's exact test for contingency table | a, b | | c, d |
The result (h,p,z) : h is true if the test rejects the null hypothesis at the alpha significance level, and false otherwise. p is the p-value and z is prior odds ratio.
runs_test ~alpha ~v x returns a test decision for the null hypothesis that the data x comes in random order, against the alternative that they do not, by running Wald–Wolfowitz runs test. The test is based on the number of runs of consecutive values above or below the mean of x. ~v is the reference value, the default value is the median of x.
mannwhitneyu ~alpha ~side x y Computes the Mann-Whitney rank test on samples x and y. If length of each sample less than 10 and no ties, then using exact test (see paper Ying Kuen Cheung and Jerome H. Klotz (1997) The Mann Whitney Wilcoxon distribution using linked list Statistica Sinica 7 805-813), else usning asymptotic normal distribution.
wilcoxon ?alpha ?side x y performs the Wilcoxon signed-rank test on the paired samples x and y.
The Wilcoxon signed-rank test is a non-parametric statistical hypothesis test used to compare two related samples, matched samples, or repeated measurements on a single sample to assess whether their population mean ranks differ.
parameteralpha
The significance level for the test, typically set to 0.05 by default. This parameter is optional.
parameterside
Specifies the alternative hypothesis. It determines whether the test is one-sided or two-sided. The default is usually a two-sided test.
parameterx
The first array of sample data.
parametery
The second array of sample data.
returns
A hypothesis type, indicating whether to reject the null hypothesis or not based on the test.
Discrete random variables
The _rvs functions generate random numbers according to the specified distribution. _pdf are "density" functions that return the probability of the element specified by the arguments, while _cdf functions are cumulative distribution functions that return the probability of all elements less than or equal to the chosen element, and _sf functions are survival functions returning one minus the corresponding CDF function. `log` versions of functions return the result for the natural logarithm of a chosen element.
binomial_rvs p n returns a random integer representing the number of successes in n trials with probability of success p on each trial.
Sourceval binomial_pdf : int ->p:float ->n:int -> float
binomial_pdf k ~p ~n returns the binomially distributed probability of k successes in n trials with probability p of success on each trial.
Sourceval binomial_logpdf : int ->p:float ->n:int -> float
binomial_logpdf k ~p ~n returns the log-binomially distributed probability of k successes in n trials with probability p of success on each trial.
Sourceval binomial_cdf : int ->p:float ->n:int -> float
binomial_cdf k ~p ~n returns the binomially distributed cumulative probability of less than or equal to k successes in n trials, with probability p on each trial.
Sourceval binomial_logcdf : int ->p:float ->n:int -> float
binomial_logcdf k ~p ~n returns the log-binomially distributed cumulative probability of less than or equal to k successes in n trials, with probability p on each trial.
Sourceval binomial_sf : int ->p:float ->n:int -> float
binomial_sf k ~p ~n is the binomial survival function, i.e. 1 - (binomial_cdf k ~p ~n).
Sourceval binomial_logsf : int ->p:float ->n:int -> float
binomial_loggf k ~p ~n is the logbinomial survival function, i.e. 1 - (binomial_logcdf k ~p ~n).
Sourceval hypergeometric_rvs : good:int ->bad:int ->sample:int -> int
hypergeometric_rvs ~good ~bad ~sample returns a random hypergeometrically distributed integer representing the number of successes in a sample (without replacement) of size ~sample from a population with ~good successful elements and ~bad unsuccessful elements.
Sourceval hypergeometric_pdf : int ->good:int ->bad:int ->sample:int -> float
hypergeometric_pdf k ~good ~bad ~sample returns the hypergeometrically distributed probability of k successes in a sample (without replacement) of ~sample elements from a population containing ~good successful elements and ~bad unsuccessful ones.
Sourceval hypergeometric_logpdf : int ->good:int ->bad:int ->sample:int -> float
hypergeometric_logpdf k ~good ~bad ~sample returns a value equivalent to a log-transformed result from hypergeometric_pdf.
Sourceval multinomial_rvs : int ->p:float array->int array
multinomial_rvs n ~p generates random numbers of multinomial distribution from n trials. The probability mass function is as follows.
p is the probability mass of k categories. The elements in p should all be positive (result is undefined if there are negative values), but they don't need to sum to 1: the result is the same whether or not p is normalized. For implementation details, refer to :cite:`davis1993computer`.
Sourceval multinomial_pdf : int array->p:float array-> float
multinomial_rvs x ~p return the probability of x given the probability mass of a multinomial distribution.
Sourceval multinomial_logpdf : int array->p:float array-> float
multinomial_rvs x ~p returns the logarithm probability of x given the probability mass of a multinomial distribution.
categorical_rvs p returns the value of a random variable which follows the categorical distribution. This is equavalent to only one trial from multinomial_rvs function, so it is just a simple wrapping.
uniform_logpdf x ~a ~b computes the natural logarithm of the probability density function (log-PDF) of the uniform distribution at the point x over the interval [a, b\).
uniform_logcdf x ~a ~b computes the natural logarithm of the cumulative distribution function (log-CDF) of the uniform distribution at the point x over the interval [a, b\).
uniform_ppf q ~a ~b computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the uniform distribution for a given probability q over the interval [a, b\).
parameterq
The probability for which to compute the corresponding quantile.
uniform_sf x ~a ~b computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the uniform distribution at the point x over the interval [a, b\).
uniform_logsf x ~a ~b computes the natural logarithm of the survival function (log-SF) of the uniform distribution at the point x over the interval [a, b\).
uniform_isf q ~a ~b computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the uniform distribution for a given probability q over the interval [a, b\).
parameterq
The probability for which to compute the corresponding value from the ISF.
exponential_logpdf x ~lambda computes the natural logarithm of the probability density function (log-PDF) of the exponential distribution at the point x with rate parameter lambda.
parameterx
The point at which to evaluate the log-PDF.
parameterlambda
The rate parameter of the exponential distribution (must be positive).
exponential_cdf x ~lambda computes the cumulative distribution function (CDF) of the exponential distribution at the point x with rate parameter lambda.
parameterx
The point at which to evaluate the CDF.
parameterlambda
The rate parameter of the exponential distribution (must be positive).
exponential_logcdf x ~lambda computes the natural logarithm of the cumulative distribution function (log-CDF) of the exponential distribution at the point x with rate parameter lambda.
parameterx
The point at which to evaluate the log-CDF.
parameterlambda
The rate parameter of the exponential distribution (must be positive).
exponential_ppf q ~lambda computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the exponential distribution for a given probability q with rate parameter lambda.
parameterq
The probability for which to compute the corresponding quantile.
parameterlambda
The rate parameter of the exponential distribution (must be positive).
exponential_sf x ~lambda computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the exponential distribution at the point x with rate parameter lambda.
parameterx
The point at which to evaluate the SF.
parameterlambda
The rate parameter of the exponential distribution (must be positive).
exponential_logsf x ~lambda computes the natural logarithm of the survival function (log-SF) of the exponential distribution at the point x with rate parameter lambda.
parameterx
The point at which to evaluate the log-SF.
parameterlambda
The rate parameter of the exponential distribution (must be positive).
exponential_isf q ~lambda computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the exponential distribution for a given probability q with rate parameter lambda.
parameterq
The probability for which to compute the corresponding value from the ISF.
parameterlambda
The rate parameter of the exponential distribution (must be positive).
exponpow_pdf x ~a ~b computes the probability density function (PDF) of the exponential power distribution at the point x with shape parameters a and b.
exponpow_logpdf x ~a ~b computes the natural logarithm of the probability density function (log-PDF) of the exponential power distribution at the point x with shape parameters a and b.
exponpow_cdf x ~a ~b computes the cumulative distribution function (CDF) of the exponential power distribution at the point x with shape parameters a and b.
exponpow_logcdf x ~a ~b computes the natural logarithm of the cumulative distribution function (log-CDF) of the exponential power distribution at the point x with shape parameters a and b.
exponpow_sf x ~a ~b computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the exponential power distribution at the point x with shape parameters a and b.
exponpow_logsf x ~a ~b computes the natural logarithm of the survival function (log-SF) of the exponential power distribution at the point x with shape parameters a and b.
gaussian_pdf x ~mu ~sigma computes the probability density function (PDF) of the Gaussian (normal) distribution at the point x with mean mu and standard deviation sigma.
gaussian_logpdf x ~mu ~sigma computes the natural logarithm of the probability density function (log-PDF) of the Gaussian (normal) distribution at the point x with mean mu and standard deviation sigma.
gaussian_cdf x ~mu ~sigma computes the cumulative distribution function (CDF) of the Gaussian (normal) distribution at the point x with mean mu and standard deviation sigma.
gaussian_logcdf x ~mu ~sigma computes the natural logarithm of the cumulative distribution function (log-CDF) of the Gaussian (normal) distribution at the point x with mean mu and standard deviation sigma.
gaussian_ppf q ~mu ~sigma computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the Gaussian (normal) distribution for a given probability q with mean mu and standard deviation sigma.
parameterq
The probability for which to compute the corresponding quantile.
gaussian_sf x ~mu ~sigma computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the Gaussian (normal) distribution at the point x with mean mu and standard deviation sigma.
gaussian_logsf x ~mu ~sigma computes the natural logarithm of the survival function (log-SF) of the Gaussian (normal) distribution at the point x with mean mu and standard deviation sigma.
gaussian_isf q ~mu ~sigma computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the Gaussian (normal) distribution for a given probability q with mean mu and standard deviation sigma.
parameterq
The probability for which to compute the corresponding value from the ISF.
gamma_pdf x ~shape ~scale computes the probability density function (PDF) of the gamma distribution at the point x with shape parameter shape and scale parameter scale.
gamma_logpdf x ~shape ~scale computes the natural logarithm of the probability density function (log-PDF) of the gamma distribution at the point x with shape parameter shape and scale parameter scale.
gamma_cdf x ~shape ~scale computes the cumulative distribution function (CDF) of the gamma distribution at the point x with shape parameter shape and scale parameter scale.
gamma_logcdf x ~shape ~scale computes the natural logarithm of the cumulative distribution function (log-CDF) of the gamma distribution at the point x with shape parameter shape and scale parameter scale.
gamma_ppf q ~shape ~scale computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the gamma distribution for a given probability q with shape parameter shape and scale parameter scale.
parameterq
The probability for which to compute the corresponding quantile.
gamma_sf x ~shape ~scale computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the gamma distribution at the point x with shape parameter shape and scale parameter scale.
gamma_logsf x ~shape ~scale computes the natural logarithm of the survival function (log-SF) of the gamma distribution at the point x with shape parameter shape and scale parameter scale.
gamma_isf q ~shape ~scale computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the gamma distribution for a given probability q with shape parameter shape and scale parameter scale.
parameterq
The probability for which to compute the corresponding value from the ISF.
beta_logpdf x ~a ~b computes the natural logarithm of the probability density function (log-PDF) of the beta distribution at the point x with shape parameters a and b.
beta_logcdf x ~a ~b computes the natural logarithm of the cumulative distribution function (log-CDF) of the beta distribution at the point x with shape parameters a and b.
beta_ppf q ~a ~b computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the beta distribution for a given probability q with shape parameters a and b.
parameterq
The probability for which to compute the corresponding quantile.
beta_sf x ~a ~b computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the beta distribution at the point x with shape parameters a and b.
beta_logsf x ~a ~b computes the natural logarithm of the survival function (log-SF) of the beta distribution at the point x with shape parameters a and b.
beta_isf q ~a ~b computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the beta distribution for a given probability q with shape parameters a and b.
parameterq
The probability for which to compute the corresponding value from the ISF.
chi2_logpdf x ~df computes the natural logarithm of the probability density function (log-PDF) of the chi-square distribution at the point x with degrees of freedom df.
chi2_logcdf x ~df computes the natural logarithm of the cumulative distribution function (log-CDF) of the chi-square distribution at the point x with degrees of freedom df.
chi2_ppf q ~df computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the chi-square distribution for a given probability q with degrees of freedom df.
parameterq
The probability for which to compute the corresponding quantile.
chi2_sf x ~df computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the chi-square distribution at the point x with degrees of freedom df.
chi2_logsf x ~df computes the natural logarithm of the survival function (log-SF) of the chi-square distribution at the point x with degrees of freedom df.
chi2_isf q ~df computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the chi-square distribution for a given probability q with degrees of freedom df.
parameterq
The probability for which to compute the corresponding value from the ISF.
f_rvs ~dfnum ~dfden generates a random variate from the F-distribution with numerator degrees of freedom dfnum and denominator degrees of freedom dfden.
The F-distribution is commonly used in the analysis of variance (ANOVA) and in the comparison of two variances.
parameterdfnum
The numerator degrees of freedom.
parameterdfden
The denominator degrees of freedom.
returns
A float representing a random sample from the F-distribution.
f_pdf x ~dfnum ~dfden computes the probability density function (PDF) of the F-distribution at the point x with numerator degrees of freedom dfnum and denominator degrees of freedom dfden.
f_logpdf x ~dfnum ~dfden computes the natural logarithm of the probability density function (log-PDF) of the F-distribution at the point x with numerator degrees of freedom dfnum and denominator degrees of freedom dfden.
f_cdf x ~dfnum ~dfden computes the cumulative distribution function (CDF) of the F-distribution at the point x with numerator degrees of freedom dfnum and denominator degrees of freedom dfden.
f_logcdf x ~dfnum ~dfden computes the natural logarithm of the cumulative distribution function (log-CDF) of the F-distribution at the point x with numerator degrees of freedom dfnum and denominator degrees of freedom dfden.
f_ppf q ~dfnum ~dfden computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the F-distribution for a given probability q with numerator degrees of freedom dfnum and denominator degrees of freedom dfden.
parameterq
The probability for which to compute the corresponding quantile.
f_sf x ~dfnum ~dfden computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the F-distribution at the point x with numerator degrees of freedom dfnum and denominator degrees of freedom dfden.
f_logsf x ~dfnum ~dfden computes the natural logarithm of the survival function (log-SF) of the F-distribution at the point x with numerator degrees of freedom dfnum and denominator degrees of freedom dfden.
f_isf q ~dfnum ~dfden computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the F-distribution for a given probability q with numerator degrees of freedom dfnum and denominator degrees of freedom dfden.
parameterq
The probability for which to compute the corresponding value from the ISF.
cauchy_pdf x ~loc ~scale computes the probability density function (PDF) of the Cauchy distribution at the point x with location parameter loc and scale parameter scale.
cauchy_logpdf x ~loc ~scale computes the natural logarithm of the probability density function (log-PDF) of the Cauchy distribution at the point x with location parameter loc and scale parameter scale.
cauchy_cdf x ~loc ~scale computes the cumulative distribution function (CDF) of the Cauchy distribution at the point x with location parameter loc and scale parameter scale.
cauchy_logcdf x ~loc ~scale computes the natural logarithm of the cumulative distribution function (log-CDF) of the Cauchy distribution at the point x with location parameter loc and scale parameter scale.
cauchy_ppf q ~loc ~scale computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the Cauchy distribution for a given probability q with location parameter loc and scale parameter scale.
parameterq
The probability for which to compute the corresponding quantile.
cauchy_sf x ~loc ~scale computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the Cauchy distribution at the point x with location parameter loc and scale parameter scale.
cauchy_logsf x ~loc ~scale computes the natural logarithm of the survival function (log-SF) of the Cauchy distribution at the point x with location parameter loc and scale parameter scale.
cauchy_isf q ~loc ~scale computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the Cauchy distribution for a given probability q with location parameter loc and scale parameter scale.
parameterq
The probability for which to compute the corresponding value from the ISF.
t_rvs ~df ~loc ~scale generates a random variate from the Student's t-distribution with df degrees of freedom, location parameter loc, and scale parameter scale.
The Student's t-distribution is commonly used in statistics for small sample sizes or when the population standard deviation is unknown.
parameterdf
The degrees of freedom of the distribution.
parameterloc
The location parameter of the distribution (the mean).
parameterscale
The scale parameter of the distribution (the standard deviation).
returns
A float representing a random sample from the t-distribution.
t_pdf x ~df ~loc ~scale computes the probability density function (PDF) of the Student's t-distribution at the point x with df degrees of freedom, location parameter loc, and scale parameter scale.
t_logpdf x ~df ~loc ~scale computes the natural logarithm of the probability density function (log-PDF) of the Student's t-distribution at the point x with df degrees of freedom, location parameter loc, and scale parameter scale.
t_cdf x ~df ~loc ~scale computes the cumulative distribution function (CDF) of the Student's t-distribution at the point x with df degrees of freedom, location parameter loc, and scale parameter scale.
t_logcdf x ~df ~loc ~scale computes the natural logarithm of the cumulative distribution function (log-CDF) of the Student's t-distribution at the point x with df degrees of freedom, location parameter loc, and scale parameter scale.
t_ppf q ~df ~loc ~scale computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the Student's t-distribution for a given probability q with df degrees of freedom, location parameter loc, and scale parameter scale.
parameterq
The probability for which to compute the corresponding quantile.
t_sf x ~df ~loc ~scale computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the Student's t-distribution at the point x with df degrees of freedom, location parameter loc, and scale parameter scale.
t_logsf x ~df ~loc ~scale computes the natural logarithm of the survival function (log-SF) of the Student's t-distribution at the point x with df degrees of freedom, location parameter loc, and scale parameter scale.
t_isf q ~df ~loc ~scale computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the Student's t-distribution for a given probability q with df degrees of freedom, location parameter loc, and scale parameter scale.
parameterq
The probability for which to compute the corresponding value from the ISF.
vonmises_pdf x ~mu ~kappa computes the probability density function (PDF) of the von Mises distribution at the point x with mean direction mu and concentration parameter kappa.
vonmises_logpdf x ~mu ~kappa computes the natural logarithm of the probability density function (log-PDF) of the von Mises distribution at the point x with mean direction mu and concentration parameter kappa.
vonmises_cdf x ~mu ~kappa computes the cumulative distribution function (CDF) of the von Mises distribution at the point x with mean direction mu and concentration parameter kappa.
vonmises_logcdf x ~mu ~kappa computes the natural logarithm of the cumulative distribution function (log-CDF) of the von Mises distribution at the point x with mean direction mu and concentration parameter kappa.
vonmises_sf x ~mu ~kappa computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the von Mises distribution at the point x with mean direction mu and concentration parameter kappa.
vonmises_logsf x ~mu ~kappa computes the natural logarithm of the survival function (log-SF) of the von Mises distribution at the point x with mean direction mu and concentration parameter kappa.
lomax_rvs ~shape ~scale generates a random variate from the Lomax distribution, also known as the Pareto distribution of the second kind, with shape parameter shape and scale parameter scale.
The Lomax distribution is often used in survival analysis and heavy-tailed modeling.
parametershape
The shape parameter of the distribution.
parameterscale
The scale parameter of the distribution.
returns
A float representing a random sample from the Lomax distribution.
lomax_pdf x ~shape ~scale computes the probability density function (PDF) of the Lomax distribution at the point x with shape parameter shape and scale parameter scale.
lomax_logpdf x ~shape ~scale computes the natural logarithm of the probability density function (log-PDF) of the Lomax distribution at the point x with shape parameter shape and scale parameter scale.
lomax_cdf x ~shape ~scale computes the cumulative distribution function (CDF) of the Lomax distribution at the point x with shape parameter shape and scale parameter scale.
lomax_logcdf x ~shape ~scale computes the natural logarithm of the cumulative distribution function (log-CDF) of the Lomax distribution at the point x with shape parameter shape and scale parameter scale.
lomax_ppf q ~shape ~scale computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the Lomax distribution for a given probability q with shape parameter shape and scale parameter scale.
parameterq
The probability for which to compute the corresponding quantile.
lomax_sf x ~shape ~scale computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the Lomax distribution at the point x with shape parameter shape and scale parameter scale.
lomax_logsf x ~shape ~scale computes the natural logarithm of the survival function (log-SF) of the Lomax distribution at the point x with shape parameter shape and scale parameter scale.
lomax_isf q ~shape ~scale computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the Lomax distribution for a given probability q with shape parameter shape and scale parameter scale.
parameterq
The probability for which to compute the corresponding value from the ISF.
weibull_pdf x ~shape ~scale computes the probability density function (PDF) of the Weibull distribution at the point x with shape parameter shape and scale parameter scale.
weibull_logpdf x ~shape ~scale computes the natural logarithm of the probability density function (log-PDF) of the Weibull distribution at the point x with shape parameter shape and scale parameter scale.
weibull_cdf x ~shape ~scale computes the cumulative distribution function (CDF) of the Weibull distribution at the point x with shape parameter shape and scale parameter scale.
weibull_logcdf x ~shape ~scale computes the natural logarithm of the cumulative distribution function (log-CDF) of the Weibull distribution at the point x with shape parameter shape and scale parameter scale.
weibull_ppf q ~shape ~scale computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the Weibull distribution for a given probability q with shape parameter shape and scale parameter scale.
parameterq
The probability for which to compute the corresponding quantile.
weibull_sf x ~shape ~scale computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the Weibull distribution at the point x with shape parameter shape and scale parameter scale.
weibull_logsf x ~shape ~scale computes the natural logarithm of the survival function (log-SF) of the Weibull distribution at the point x with shape parameter shape and scale parameter scale.
weibull_isf q ~shape ~scale computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the Weibull distribution for a given probability q with shape parameter shape and scale parameter scale.
parameterq
The probability for which to compute the corresponding value from the ISF.
laplace_rvs ~loc ~scale generates a random variate from the Laplace distribution, also known as the double exponential distribution, with location parameter loc and scale parameter scale.
The Laplace distribution is often used in statistical methods for detecting outliers.
parameterloc
The location parameter of the distribution (the mean).
parameterscale
The scale parameter of the distribution (the standard deviation).
returns
A float representing a random sample from the Laplace distribution.
laplace_pdf x ~loc ~scale computes the probability density function (PDF) of the Laplace distribution at the point x with location parameter loc and scale parameter scale.
laplace_logpdf x ~loc ~scale computes the natural logarithm of the probability density function (log-PDF) of the Laplace distribution at the point x with location parameter loc and scale parameter scale.
laplace_cdf x ~loc ~scale computes the cumulative distribution function (CDF) of the Laplace distribution at the point x with location parameter loc and scale parameter scale.
laplace_logcdf x ~loc ~scale computes the natural logarithm of the cumulative distribution function (log-CDF) of the Laplace distribution at the point x with location parameter loc and scale parameter scale.
laplace_ppf q ~loc ~scale computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the Laplace distribution for a given probability q with location parameter loc and scale parameter scale.
parameterq
The probability for which to compute the corresponding quantile.
laplace_sf x ~loc ~scale computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the Laplace distribution at the point x with location parameter loc and scale parameter scale.
laplace_logsf x ~loc ~scale computes the natural logarithm of the survival function (log-SF) of the Laplace distribution at the point x with location parameter loc and scale parameter scale.
laplace_isf q ~loc ~scale computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the Laplace distribution for a given probability q with location parameter loc and scale parameter scale.
parameterq
The probability for which to compute the corresponding value from the ISF.
gumbel1_pdf x ~a ~b computes the probability density function (PDF) of the Gumbel Type I distribution at the point x with location parameter a and scale parameter b.
gumbel1_logpdf x ~a ~b computes the natural logarithm of the probability density function (log-PDF) of the Gumbel Type I distribution at the point x with location parameter a and scale parameter b.
gumbel1_cdf x ~a ~b computes the cumulative distribution function (CDF) of the Gumbel Type I distribution at the point x with location parameter a and scale parameter b.
gumbel1_logcdf x ~a ~b computes the natural logarithm of the cumulative distribution function (log-CDF) of the Gumbel Type I distribution at the point x with location parameter a and scale parameter b.
gumbel1_ppf q ~a ~b computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the Gumbel Type I distribution for a given probability q with location parameter a and scale parameter b.
parameterq
The probability for which to compute the corresponding quantile.
gumbel1_sf x ~a ~b computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the Gumbel Type I distribution at the point x with location parameter a and scale parameter b.
gumbel1_logsf x ~a ~b computes the natural logarithm of the survival function (log-SF) of the Gumbel Type I distribution at the point x with location parameter a and scale parameter b.
gumbel1_isf q ~a ~b computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the Gumbel Type I distribution for a given probability q with location parameter a and scale parameter b.
parameterq
The probability for which to compute the corresponding value from the ISF.
gumbel2_pdf x ~a ~b computes the probability density function (PDF) of the Gumbel Type II distribution at the point x with location parameter a and scale parameter b.
gumbel2_logpdf x ~a ~b computes the natural logarithm of the probability density function (log-PDF) of the Gumbel Type II distribution at the point x with location parameter a and scale parameter b.
gumbel2_cdf x ~a ~b computes the cumulative distribution function (CDF) of the Gumbel Type II distribution at the point x with location parameter a and scale parameter b.
gumbel2_logcdf x ~a ~b computes the natural logarithm of the cumulative distribution function (log-CDF) of the Gumbel Type II distribution at the point x with location parameter a and scale parameter b.
gumbel2_ppf q ~a ~b computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the Gumbel Type II distribution for a given probability q with location parameter a and scale parameter b.
parameterq
The probability for which to compute the corresponding quantile.
gumbel2_sf x ~a ~b computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the Gumbel Type II distribution at the point x with location parameter a and scale parameter b.
gumbel2_logsf x ~a ~b computes the natural logarithm of the survival function (log-SF) of the Gumbel Type II distribution at the point x with location parameter a and scale parameter b.
gumbel2_isf q ~a ~b computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the Gumbel Type II distribution for a given probability q with location parameter a and scale parameter b.
parameterq
The probability for which to compute the corresponding value from the ISF.
logistic_pdf x ~loc ~scale computes the probability density function (PDF) of the logistic distribution at the point x with location parameter loc and scale parameter scale.
logistic_logpdf x ~loc ~scale computes the natural logarithm of the probability density function (log-PDF) of the logistic distribution at the point x with location parameter loc and scale parameter scale.
logistic_cdf x ~loc ~scale computes the cumulative distribution function (CDF) of the logistic distribution at the point x with location parameter loc and scale parameter scale.
logistic_logcdf x ~loc ~scale computes the natural logarithm of the cumulative distribution function (log-CDF) of the logistic distribution at the point x with location parameter loc and scale parameter scale.
logistic_ppf q ~loc ~scale computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the logistic distribution for a given probability q with location parameter loc and scale parameter scale.
parameterq
The probability for which to compute the corresponding quantile.
logistic_sf x ~loc ~scale computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the logistic distribution at the point x with location parameter loc and scale parameter scale.
logistic_logsf x ~loc ~scale computes the natural logarithm of the survival function (log-SF) of the logistic distribution at the point x with location parameter loc and scale parameter scale.
logistic_isf q ~loc ~scale computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the logistic distribution for a given probability q with location parameter loc and scale parameter scale.
parameterq
The probability for which to compute the corresponding value from the ISF.
lognormal_rvs ~mu ~sigma generates a random variate from the log-normal distribution with parameters mu (mean of the underlying normal distribution) and sigma (standard deviation of the underlying normal distribution).
The log-normal distribution is commonly used to model positive-valued data with a distribution that is skewed to the right.
parametermu
The mean of the underlying normal distribution.
parametersigma
The standard deviation of the underlying normal distribution.
returns
A float representing a random sample from the log-normal distribution.
lognormal_pdf x ~mu ~sigma computes the probability density function (PDF) of the log-normal distribution at the point x with parameters mu (mean of the underlying normal distribution) and sigma (standard deviation of the underlying normal distribution).
parameterx
The point at which to evaluate the PDF.
parametermu
The mean of the underlying normal distribution.
parametersigma
The standard deviation of the underlying normal distribution.
lognormal_logpdf x ~mu ~sigma computes the natural logarithm of the probability density function (log-PDF) of the log-normal distribution at the point x with parameters mu and sigma.
parameterx
The point at which to evaluate the log-PDF.
parametermu
The mean of the underlying normal distribution.
parametersigma
The standard deviation of the underlying normal distribution.
lognormal_cdf x ~mu ~sigma computes the cumulative distribution function (CDF) of the log-normal distribution at the point x with parameters mu and sigma.
parameterx
The point at which to evaluate the CDF.
parametermu
The mean of the underlying normal distribution.
parametersigma
The standard deviation of the underlying normal distribution.
lognormal_logcdf x ~mu ~sigma computes the natural logarithm of the cumulative distribution function (log-CDF) of the log-normal distribution at the point x with parameters mu and sigma.
parameterx
The point at which to evaluate the log-CDF.
parametermu
The mean of the underlying normal distribution.
parametersigma
The standard deviation of the underlying normal distribution.
lognormal_ppf q ~mu ~sigma computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the log-normal distribution for a given probability q with parameters mu and sigma.
parameterq
The probability for which to compute the corresponding quantile.
parametermu
The mean of the underlying normal distribution.
parametersigma
The standard deviation of the underlying normal distribution.
lognormal_sf x ~mu ~sigma computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the log-normal distribution at the point x with parameters mu and sigma.
parameterx
The point at which to evaluate the SF.
parametermu
The mean of the underlying normal distribution.
parametersigma
The standard deviation of the underlying normal distribution.
lognormal_logsf x ~mu ~sigma computes the natural logarithm of the survival function (log-SF) of the log-normal distribution at the point x with parameters mu and sigma.
parameterx
The point at which to evaluate the log-SF.
parametermu
The mean of the underlying normal distribution.
parametersigma
The standard deviation of the underlying normal distribution.
lognormal_isf q ~mu ~sigma computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the log-normal distribution for a given probability q with parameters mu and sigma.
parameterq
The probability for which to compute the corresponding value from the ISF.
parametermu
The mean of the underlying normal distribution.
parametersigma
The standard deviation of the underlying normal distribution.
rayleigh_logpdf x ~sigma computes the natural logarithm of the probability density function (log-PDF) of the Rayleigh distribution at the point x with scale parameter sigma.
rayleigh_logcdf x ~sigma computes the natural logarithm of the cumulative distribution function (log-CDF) of the Rayleigh distribution at the point x with scale parameter sigma.
rayleigh_ppf q ~sigma computes the percent-point function (PPF), also known as the quantile function or inverse CDF, of the Rayleigh distribution for a given probability q with scale parameter sigma.
parameterq
The probability for which to compute the corresponding quantile.
rayleigh_sf x ~sigma computes the survival function (SF), which is one minus the cumulative distribution function (1 - CDF), of the Rayleigh distribution at the point x with scale parameter sigma.
rayleigh_logsf x ~sigma computes the natural logarithm of the survival function (log-SF) of the Rayleigh distribution at the point x with scale parameter sigma.
rayleigh_isf q ~sigma computes the inverse survival function (ISF), which is the inverse of the survival function (SF), of the Rayleigh distribution for a given probability q with scale parameter sigma.
parameterq
The probability for which to compute the corresponding value from the ISF.
dirichlet_pdf x ~alpha computes the probability density function (PDF) of the Dirichlet distribution for the input vector x with concentration parameters alpha.
The Dirichlet distribution is a multivariate generalization of the Beta distribution and is commonly used as a prior distribution in Bayesian statistics, particularly in the context of categorical data and multinomial distributions.
parameterx
The input vector for which to evaluate the PDF, typically representing proportions that sum to 1.
parameteralpha
The concentration parameters of the distribution, where each element must be greater than 0.
dirichlet_logpdf x ~alpha computes the natural logarithm of the probability density function (log-PDF) of the Dirichlet distribution for the input vector x with concentration parameters alpha.
parameterx
The input vector for which to evaluate the log-PDF, typically representing proportions that sum to 1.
parameteralpha
The concentration parameters of the distribution, where each element must be greater than 0.