package owl

  1. Overview
  2. Docs
Legend:
Page
Library
Module
Module type
Parameter
Class
Class type
Source

Module Owl_statsSource

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.

Randomisation functions
Sourceval shuffle : 'a array -> 'a array

``shuffle x`` return a new array of the shuffled ``x``.

Sourceval choose : 'a array -> int -> 'a array

``choose x n`` draw ``n`` samples from ``x`` without replecement.

Sourceval sample : 'a array -> int -> 'a array

``sample x n`` draw ``n`` samples from ``x`` with replacement.

Basic statistical functions
Sourceval sum : float array -> float

``sum x`` returns the summation of the elements in ``x``.

Sourceval mean : float array -> float

``mean x`` returns the mean of the elements in ``x``.

Sourceval var : ?mean:float -> float array -> float

``var x`` returns the variance of elements in ``x``.

Sourceval std : ?mean:float -> float array -> float

``std x`` calculates the standard deviation of ``x``.

Sourceval sem : ?mean:float -> float array -> float

``sem x`` calculates the standard error of ``x``, also referred to as standard error of the mean.

Sourceval absdev : ?mean:float -> float array -> float

``absdev x`` calculates the average absolute deviation of ``x``.

Sourceval skew : ?mean:float -> ?sd:float -> float array -> float

``skew x`` calculates the skewness (the third standardized moment) of ``x``.

Sourceval kurtosis : ?mean:float -> ?sd:float -> float array -> float

``kurtosis x`` calculates the Pearson's kurtosis of ``x``, i.e. the fourth standardized moment of ``x``.

Sourceval central_moment : int -> float array -> float

``central_moment n x`` calculates the ``n`` th central moment of ``x``.

Sourceval cov : ?m0:float -> ?m1:float -> float array -> float array -> float

``cov x0 x1`` calculates the covariance of ``x0`` and ``x1``, the mean of ``x0`` and ``x1`` can be specified by ``m0`` and ``m1`` respectively.

Sourceval concordant : 'a array -> 'b array -> int

TODO

Sourceval discordant : 'a array -> 'b array -> int

TODO

Sourceval corrcoef : float array -> float array -> float

``corrcoef x y`` calculates the Pearson correlation of ``x`` and ``y``. Namely, ``corrcoef x y = cov(x, y) / (sigma_x * sigma_y)``.

Sourceval kendall_tau : float array -> float array -> float

``kendall_tau x y`` calculates the Kendall Tau correlation between ``x`` and ``y``.

Sourceval spearman_rho : float array -> float array -> float

``spearman_rho x y`` calculates the Spearman Rho correlation between ``x`` and ``y``.

Sourceval autocorrelation : ?lag:int -> float array -> float

``autocorrelation ~lag x`` calculates the autocorrelation of ``x`` with the given ``lag``.

Sourceval percentile : float array -> float -> float

``percentile x p`` returns the ``p`` percentile of the data ``x``. ``p`` is between 0. and 100. ``x`` does not need to be sorted beforehand.

Sourceval quantile : float array -> float -> float

``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.

Sourceval first_quartile : float array -> float

``first_quartile x`` returns the first quartile of ``x``, i.e. 25 percentiles.

Sourceval third_quartile : float array -> float

``third_quartile x`` returns the third quartile of ``x``, i.e. 75 percentiles.

Sourceval interquartile : float array -> float

``interquartile x`` returns the interquartile range of ``x`` which is defined as the first quartile subtracted from the third quartile.

Sourceval median : float array -> float

``median x`` returns the median of ``x``.

Sourceval min : float array -> float

``min x`` returns the minimum element in ``x``.

Sourceval max : float array -> float

``max x`` returns the maximum element in ``x``.

Sourceval minmax : float array -> float * float

``minmax x`` returns both ``(minimum, maximum)`` elements in ``x``.

Sourceval min_i : float array -> int

``min_i x`` returns the index of the minimum in ``x``.

Sourceval max_i : float array -> int

``max_i x`` returns the index of the maximum in ``x``.

Sourceval minmax_i : float array -> int * int

``minmax_i x`` returns the indices of both minimum and maximum in ``x``.

Sourceval sort : ?inc:bool -> float array -> float array

``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.

Sourceval argsort : ?inc:bool -> float array -> int array

``argsort x`` sorts the elements in ``x`` and returns the indices mapping of the elements in the current array to their original position in ``x``.

The sorting is in increasing order if ``inc = true``, otherwise in decreasing order if ``inc=false``. By default, ``inc`` is ``true``.

Sourceval rank : ?ties_strategy:[ `Average | `Min | `Max ] -> float array -> float array

Computes sample's ranks.

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.

Sourceval normalise : histogram -> histogram

``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.

Sourceval normalise_density : histogram -> histogram

``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.

Sourceval pp_hist : Format.formatter -> histogram -> unit

Pretty-print summary information on a histogram record

Sourceval ecdf : float array -> float array * float array

``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.

Sourceval z_score : mu:float -> sigma:float -> float array -> float array

``z_score x`` calculates the z score of a given array ``x``.

Sourceval t_score : float array -> float array

``t_score x`` calculates the t score of a given array ``x``.

Sourceval normlise_pdf : float array -> float array

TODO

Sourceval tukey_fences : ?k:float -> float array -> float * float

``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:

.. math:: (Q1 - k*(Q3-Q1), Q3 + k*(Q3-Q1))

Sourceval gaussian_kde : ?bandwidth:[ `Silverman | `Scott ] -> ?n_points:int -> float array -> float array * float array

``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.

The default bandwidth value is ``Scott``.

MCMC: Markov Chain Monte Carlo
Sourceval metropolis_hastings : (float array -> float) -> float array -> int -> float array array

TODO: ``metropolis_hastings f p n`` is Metropolis-Hastings MCMC algorithm. f is pdf of the p

Sourceval gibbs_sampling : (float array -> int -> float) -> float array -> int -> float array array

TODO: ``gibbs_sampling f p n`` is Gibbs sampler. f is a sampler based on the full conditional function of all variables

Hypothesis tests
Sourcetype hypothesis = {
  1. reject : bool;
  2. p_value : float;
  3. score : float;
}

Record type contains the result of a hypothesis test.

Sourcetype tail =
  1. | BothSide
  2. | RightSide
  3. | LeftSide
    (*

    Types of alternative hypothesis tests: one-side, left-side, or right-side.

    *)
Sourceval z_test : mu:float -> sigma:float -> ?alpha:float -> ?side:tail -> float array -> hypothesis

``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.

Sourceval t_test : mu:float -> ?alpha:float -> ?side:tail -> float array -> hypothesis

``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.

Sourceval t_test_paired : ?alpha:float -> ?side:tail -> float array -> float array -> hypothesis

``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.

Sourceval t_test_unpaired : ?alpha:float -> ?side:tail -> ?equal_var:bool -> float array -> float array -> hypothesis

``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.

Sourceval ks_test : ?alpha:float -> float array -> (float -> float) -> hypothesis

``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.

Sourceval ks2_test : ?alpha:float -> float array -> float array -> hypothesis

``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.

Sourceval var_test : ?alpha:float -> ?side:tail -> variance:float -> float array -> hypothesis

``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.

Sourceval jb_test : ?alpha:float -> float array -> hypothesis

``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.

Sourceval fisher_test : ?alpha:float -> ?side:tail -> int -> int -> int -> int -> hypothesis

``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.

Sourceval runs_test : ?alpha:float -> ?side:tail -> ?v:float -> float array -> hypothesis

``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 runnign 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``.

Sourceval mannwhitneyu : ?alpha:float -> ?side:tail -> float array -> float array -> hypothesis

``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.

Sourceval wilcoxon : ?alpha:float -> ?side:tail -> float array -> float array -> hypothesis

TODO

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.

Sourceval uniform_int_rvs : a:int -> b:int -> int

``uniform_rvs ~a ~b`` returns a random uniformly distributed integer between ``a`` and ``b``, inclusive.

Sourceval binomial_rvs : p:float -> n:int -> int

``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 probabilty mass function is as follows.

.. math:: P(x) = \fracn!{x_1! \cdot\cdot\cdot x_k!

}

p_

^x_1 \cdot\cdot\cdot p_k^x_k

``p`` is the probabilty mass of ``k`` categories. If the elements in ``p`` do not sum to 1, the last element of the ``p`` array is not used and is replaced with the remaining probability left over from the earlier elements.

For implemantation, refer to :cite:`davis1993computer`.

Sourceval multinomial_pdf : int array -> p:float array -> float

``multinomial_rvs x ~p`` return the probabilty of ``x`` given the probabilty mass of a multinomial distribution.

Sourceval multinomial_logpdf : int array -> p:float array -> float

``multinomial_rvs x ~p`` returns the logarithm probabilty of ``x`` given the probabilty mass of a multinomial distribution.

Sourceval categorical_rvs : float array -> int

``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.

Continuous random variables
Sourceval std_uniform_rvs : unit -> float

TODO

Sourceval uniform_rvs : a:float -> b:float -> float

TODO

Sourceval uniform_pdf : float -> a:float -> b:float -> float

TODO

Sourceval uniform_logpdf : float -> a:float -> b:float -> float

TODO

Sourceval uniform_cdf : float -> a:float -> b:float -> float

TODO

Sourceval uniform_logcdf : float -> a:float -> b:float -> float

TODO

Sourceval uniform_ppf : float -> a:float -> b:float -> float

TODO

Sourceval uniform_sf : float -> a:float -> b:float -> float

TODO

Sourceval uniform_logsf : float -> a:float -> b:float -> float

TODO

Sourceval uniform_isf : float -> a:float -> b:float -> float

TODO

Sourceval exponential_rvs : lambda:float -> float

TODO

Sourceval exponential_pdf : float -> lambda:float -> float

TODO

Sourceval exponential_logpdf : float -> lambda:float -> float

TODO

Sourceval exponential_cdf : float -> lambda:float -> float

TODO

Sourceval exponential_logcdf : float -> lambda:float -> float

TODO

Sourceval exponential_ppf : float -> lambda:float -> float

TODO

Sourceval exponential_sf : float -> lambda:float -> float

TODO

Sourceval exponential_logsf : float -> lambda:float -> float

TODO

Sourceval exponential_isf : float -> lambda:float -> float

TODO

Sourceval exponpow_rvs : a:float -> b:float -> float

.. math:: p(x) dx = (1/(2 a Gamma(1+1/b))) * exp(-|x/a|^b) dx

Sourceval exponpow_pdf : float -> a:float -> b:float -> float

TODO

Sourceval exponpow_logpdf : float -> a:float -> b:float -> float

TODO

Sourceval exponpow_cdf : float -> a:float -> b:float -> float

TODO

Sourceval exponpow_logcdf : float -> a:float -> b:float -> float

TODO

Sourceval exponpow_sf : float -> a:float -> b:float -> float

TODO

Sourceval exponpow_logsf : float -> a:float -> b:float -> float

TODO

Sourceval gaussian_rvs : mu:float -> sigma:float -> float

TODO

Sourceval gaussian_pdf : float -> mu:float -> sigma:float -> float

TODO

Sourceval gaussian_logpdf : float -> mu:float -> sigma:float -> float

TODO

Sourceval gaussian_cdf : float -> mu:float -> sigma:float -> float

TODO

Sourceval gaussian_logcdf : float -> mu:float -> sigma:float -> float

TODO

Sourceval gaussian_ppf : float -> mu:float -> sigma:float -> float

TODO

Sourceval gaussian_sf : float -> mu:float -> sigma:float -> float

TODO

Sourceval gaussian_logsf : float -> mu:float -> sigma:float -> float

TODO

Sourceval gaussian_isf : float -> mu:float -> sigma:float -> float

TODO

Sourceval gamma_rvs : shape:float -> scale:float -> float

TODO

Sourceval gamma_pdf : float -> shape:float -> scale:float -> float

TODO

Sourceval gamma_logpdf : float -> shape:float -> scale:float -> float

TODO

Sourceval gamma_cdf : float -> shape:float -> scale:float -> float

TODO

Sourceval gamma_logcdf : float -> shape:float -> scale:float -> float

TODO

Sourceval gamma_ppf : float -> shape:float -> scale:float -> float

TODO

Sourceval gamma_sf : float -> shape:float -> scale:float -> float

TODO

Sourceval gamma_logsf : float -> shape:float -> scale:float -> float

TODO

Sourceval gamma_isf : float -> shape:float -> scale:float -> float

TODO

Sourceval beta_rvs : a:float -> b:float -> float

TODO

Sourceval beta_pdf : float -> a:float -> b:float -> float

TODO

Sourceval beta_logpdf : float -> a:float -> b:float -> float

TODO

Sourceval beta_cdf : float -> a:float -> b:float -> float

TODO

Sourceval beta_logcdf : float -> a:float -> b:float -> float

TODO

Sourceval beta_ppf : float -> a:float -> b:float -> float

TODO

Sourceval beta_sf : float -> a:float -> b:float -> float

TODO

Sourceval beta_logsf : float -> a:float -> b:float -> float

TODO

Sourceval beta_isf : float -> a:float -> b:float -> float

TODO

Sourceval chi2_rvs : df:float -> float

TODO

Sourceval chi2_pdf : float -> df:float -> float

TODO

Sourceval chi2_logpdf : float -> df:float -> float

TODO

Sourceval chi2_cdf : float -> df:float -> float

TODO

Sourceval chi2_logcdf : float -> df:float -> float

TODO

Sourceval chi2_ppf : float -> df:float -> float

TODO

Sourceval chi2_sf : float -> df:float -> float

TODO

Sourceval chi2_logsf : float -> df:float -> float

TODO

Sourceval chi2_isf : float -> df:float -> float

TODO

Sourceval f_rvs : dfnum:float -> dfden:float -> float

TODO

Sourceval f_pdf : float -> dfnum:float -> dfden:float -> float

TODO

Sourceval f_logpdf : float -> dfnum:float -> dfden:float -> float

TODO

Sourceval f_cdf : float -> dfnum:float -> dfden:float -> float

TODO

Sourceval f_logcdf : float -> dfnum:float -> dfden:float -> float

TODO

Sourceval f_ppf : float -> dfnum:float -> dfden:float -> float

TODO

Sourceval f_sf : float -> dfnum:float -> dfden:float -> float

TODO

Sourceval f_logsf : float -> dfnum:float -> dfden:float -> float

TODO

Sourceval f_isf : float -> dfnum:float -> dfden:float -> float

TODO

Sourceval cauchy_rvs : loc:float -> scale:float -> float

TODO

Sourceval cauchy_pdf : float -> loc:float -> scale:float -> float

TODO

Sourceval cauchy_logpdf : float -> loc:float -> scale:float -> float

TODO

Sourceval cauchy_cdf : float -> loc:float -> scale:float -> float

TODO

Sourceval cauchy_logcdf : float -> loc:float -> scale:float -> float

TODO

Sourceval cauchy_ppf : float -> loc:float -> scale:float -> float

TODO

Sourceval cauchy_sf : float -> loc:float -> scale:float -> float

TODO

Sourceval cauchy_logsf : float -> loc:float -> scale:float -> float

TODO

Sourceval cauchy_isf : float -> loc:float -> scale:float -> float

TODO

Sourceval t_rvs : df:float -> loc:float -> scale:float -> float

TODO

Sourceval t_pdf : float -> df:float -> loc:float -> scale:float -> float

TODO

Sourceval t_logpdf : float -> df:float -> loc:float -> scale:float -> float

TODO

Sourceval t_cdf : float -> df:float -> loc:float -> scale:float -> float

TODO

Sourceval t_logcdf : float -> df:float -> loc:float -> scale:float -> float

TODO

Sourceval t_ppf : float -> df:float -> loc:float -> scale:float -> float

TODO

Sourceval t_sf : float -> df:float -> loc:float -> scale:float -> float

TODO

Sourceval t_logsf : float -> df:float -> loc:float -> scale:float -> float

TODO

Sourceval t_isf : float -> df:float -> loc:float -> scale:float -> float

TODO

Sourceval vonmises_rvs : mu:float -> kappa:float -> float

TODO

Sourceval vonmises_pdf : float -> mu:float -> kappa:float -> float

TODO

Sourceval vonmises_logpdf : float -> mu:float -> kappa:float -> float

TODO

Sourceval vonmises_cdf : float -> mu:float -> kappa:float -> float

TODO

Sourceval vonmises_logcdf : float -> mu:float -> kappa:float -> float

TODO

Sourceval vonmises_sf : float -> mu:float -> kappa:float -> float

TODO

Sourceval vonmises_logsf : float -> mu:float -> kappa:float -> float

TODO

Sourceval lomax_rvs : shape:float -> scale:float -> float

TODO

Sourceval lomax_pdf : float -> shape:float -> scale:float -> float

TODO

Sourceval lomax_logpdf : float -> shape:float -> scale:float -> float

TODO

Sourceval lomax_cdf : float -> shape:float -> scale:float -> float

TODO

Sourceval lomax_logcdf : float -> shape:float -> scale:float -> float

TODO

Sourceval lomax_ppf : float -> shape:float -> scale:float -> float

TODO

Sourceval lomax_sf : float -> shape:float -> scale:float -> float

TODO

Sourceval lomax_logsf : float -> shape:float -> scale:float -> float

TODO

Sourceval lomax_isf : float -> shape:float -> scale:float -> float

TODO

Sourceval weibull_rvs : shape:float -> scale:float -> float

TODO

Sourceval weibull_pdf : float -> shape:float -> scale:float -> float

TODO

Sourceval weibull_logpdf : float -> shape:float -> scale:float -> float

TODO

Sourceval weibull_cdf : float -> shape:float -> scale:float -> float

TODO

Sourceval weibull_logcdf : float -> shape:float -> scale:float -> float

TODO

Sourceval weibull_ppf : float -> shape:float -> scale:float -> float

TODO

Sourceval weibull_sf : float -> shape:float -> scale:float -> float

TODO

Sourceval weibull_logsf : float -> shape:float -> scale:float -> float

TODO

Sourceval weibull_isf : float -> shape:float -> scale:float -> float

TODO

Sourceval laplace_rvs : loc:float -> scale:float -> float

TODO

Sourceval laplace_pdf : float -> loc:float -> scale:float -> float

TODO

Sourceval laplace_logpdf : float -> loc:float -> scale:float -> float

TODO

Sourceval laplace_cdf : float -> loc:float -> scale:float -> float

TODO

Sourceval laplace_logcdf : float -> loc:float -> scale:float -> float

TODO

Sourceval laplace_ppf : float -> loc:float -> scale:float -> float

TODO

Sourceval laplace_sf : float -> loc:float -> scale:float -> float

TODO

Sourceval laplace_logsf : float -> loc:float -> scale:float -> float

TODO

Sourceval laplace_isf : float -> loc:float -> scale:float -> float

TODO

Sourceval gumbel1_rvs : a:float -> b:float -> float

TODO

Sourceval gumbel1_pdf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel1_logpdf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel1_cdf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel1_logcdf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel1_ppf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel1_sf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel1_logsf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel1_isf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel2_rvs : a:float -> b:float -> float

TODO

Sourceval gumbel2_pdf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel2_logpdf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel2_cdf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel2_logcdf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel2_ppf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel2_sf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel2_logsf : float -> a:float -> b:float -> float

TODO

Sourceval gumbel2_isf : float -> a:float -> b:float -> float

TODO

Sourceval logistic_rvs : loc:float -> scale:float -> float

TODO

Sourceval logistic_pdf : float -> loc:float -> scale:float -> float

TODO

Sourceval logistic_logpdf : float -> loc:float -> scale:float -> float

TODO

Sourceval logistic_cdf : float -> loc:float -> scale:float -> float

TODO

Sourceval logistic_logcdf : float -> loc:float -> scale:float -> float

TODO

Sourceval logistic_ppf : float -> loc:float -> scale:float -> float

TODO

Sourceval logistic_sf : float -> loc:float -> scale:float -> float

TODO

Sourceval logistic_logsf : float -> loc:float -> scale:float -> float

TODO

Sourceval logistic_isf : float -> loc:float -> scale:float -> float

TODO

Sourceval lognormal_rvs : mu:float -> sigma:float -> float

TODO

Sourceval lognormal_pdf : float -> mu:float -> sigma:float -> float

TODO

Sourceval lognormal_logpdf : float -> mu:float -> sigma:float -> float

TODO

Sourceval lognormal_cdf : float -> mu:float -> sigma:float -> float

TODO

Sourceval lognormal_logcdf : float -> mu:float -> sigma:float -> float

TODO

Sourceval lognormal_ppf : float -> mu:float -> sigma:float -> float

TODO

Sourceval lognormal_sf : float -> mu:float -> sigma:float -> float

TODO

Sourceval lognormal_logsf : float -> mu:float -> sigma:float -> float

TODO

Sourceval lognormal_isf : float -> mu:float -> sigma:float -> float

TODO

Sourceval rayleigh_rvs : sigma:float -> float

TODO

Sourceval rayleigh_pdf : float -> sigma:float -> float

TODO

Sourceval rayleigh_logpdf : float -> sigma:float -> float

TODO

Sourceval rayleigh_cdf : float -> sigma:float -> float

TODO

Sourceval rayleigh_logcdf : float -> sigma:float -> float

TODO

Sourceval rayleigh_ppf : float -> sigma:float -> float

TODO

Sourceval rayleigh_sf : float -> sigma:float -> float

TODO

Sourceval rayleigh_logsf : float -> sigma:float -> float

TODO

Sourceval rayleigh_isf : float -> sigma:float -> float

TODO

Sourceval dirichlet_rvs : alpha:float array -> float array

``dirichlet_rvs ~alpha`` returns random variables of ``K-1`` order Dirichlet distribution, follows the following probabilty dense function.

.. math:: f(x_1,...,x_K; \alpha_1,...,\alpha_K) = \frac

\mathbf{B(\alpha)

}

\prod_=1^K x_i^\alpha_i - 1

The normalising constant is the multivariate Beta function, which can be expressed in terms of the gamma function:

.. math:: \mathbfB(\alpha) = \frac\prod_{i=1^K \Gamma(\alpha_i)

}

\Gamma(\sum_{i=1^K \alpha_i)

}

Note that ``x`` is a standard K-1 simplex, i.e. :math:`\sum_i^K x_i = 1` and :math:`x_i \ge 0, \forall x_i \in 1,K`.

Sourceval dirichlet_pdf : float array -> alpha:float array -> float

TODO

Sourceval dirichlet_logpdf : float array -> alpha:float array -> float

TODO