Class InferenceTestUtils


  • public class InferenceTestUtils
    extends Object
    A collection of static methods to create inference test instances or to perform inference tests.
    • Method Summary

      All Methods Static Methods Concrete Methods 
      Modifier and Type Method Description
      static double approximateP​(double d, int n, int m)
      Uses the Kolmogorov-Smirnov distribution to approximate \(P(D_{n,m} > d)\) where \(D_{n,m}\) is the 2-sample Kolmogorov-Smirnov statistic.
      static double chiSquare​(double[] expected, long[] observed)
      Computes the Chi-Square statistic comparing observed and expected frequency counts.
      static double chiSquare​(long[][] counts)
      Computes the Chi-Square statistic associated with a chi-square test of independence based on the input counts array, viewed as a two-way table.
      static double chiSquareDataSetsComparison​(long[] observed1, long[] observed2)
      Computes a Chi-Square two sample test statistic comparing bin frequency counts in observed1 and observed2.
      static double chiSquareTest​(double[] expected, long[] observed)
      Returns the observed significance level, or p-value, associated with a Chi-square goodness of fit test comparing the observed frequency counts to those in the expected array.
      static boolean chiSquareTest​(double[] expected, long[] observed, double alpha)
      Performs a Chi-square goodness of fit test evaluating the null hypothesis that the observed counts conform to the frequency distribution described by the expected counts, with significance level alpha.
      static double chiSquareTest​(long[][] counts)
      Returns the observed significance level, or p-value, associated with a chi-square test of independence based on the input counts array, viewed as a two-way table.
      static boolean chiSquareTest​(long[][] counts, double alpha)
      Performs a chi-square test of independence evaluating the null hypothesis that the classifications represented by the counts in the columns of the input 2-way table are independent of the rows, with significance level alpha.
      static double chiSquareTestDataSetsComparison​(long[] observed1, long[] observed2)
      Returns the observed significance level, or p-value, associated with a Chi-Square two sample test comparing bin frequency counts in observed1 and observed2.
      static boolean chiSquareTestDataSetsComparison​(long[] observed1, long[] observed2, double alpha)
      Performs a Chi-Square two sample test comparing two binned data sets.
      static double exactP​(double d, int m, int n, boolean strict)
      Computes \(P(D_{n,m} > d)\) if strict is true; otherwise \(P(D_{n,m} \ge d)\), where \(D_{n,m}\) is the 2-sample Kolmogorov-Smirnov statistic.
      static double g​(double[] expected, long[] observed)
      Computes the G statistic for Goodness of Fit comparing observed and expected frequency counts.
      static double gDataSetsComparison​(long[] observed1, long[] observed2)
      Computes a G (Log-Likelihood Ratio) two sample test statistic for independence comparing frequency counts in observed1 and observed2.
      static double gTest​(double[] expected, long[] observed)
      Returns the observed significance level, or p-value, associated with a G-Test for goodness of fit comparing the observed frequency counts to those in the expected array.
      static boolean gTest​(double[] expected, long[] observed, double alpha)
      Performs a G-Test (Log-Likelihood Ratio Test) for goodness of fit evaluating the null hypothesis that the observed counts conform to the frequency distribution described by the expected counts, with significance level alpha.
      static double gTestDataSetsComparison​(long[] observed1, long[] observed2)
      Returns the observed significance level, or p-value, associated with a G-Value (Log-Likelihood Ratio) for two sample test comparing bin frequency counts in observed1 and observed2.
      static boolean gTestDataSetsComparison​(long[] observed1, long[] observed2, double alpha)
      Performs a G-Test (Log-Likelihood Ratio Test) comparing two binned data sets.
      static double gTestIntrinsic​(double[] expected, long[] observed)
      Returns the intrinsic (Hardy-Weinberg proportions) p-Value, as described in p64-69 of McDonald, J.H. 2009.
      static double homoscedasticT​(double[] sample1, double[] sample2)
      Computes a 2-sample t statistic, under the hypothesis of equal subpopulation variances.
      static double homoscedasticT​(StatisticalSummary sampleStats1, StatisticalSummary sampleStats2)
      Computes a 2-sample t statistic, comparing the means of the datasets described by two StatisticalSummary instances, under the assumption of equal subpopulation variances.
      static double homoscedasticTTest​(double[] sample1, double[] sample2)
      Returns the observed significance level, or p-value, associated with a two-sample, two-tailed t-test comparing the means of the input arrays, under the assumption that the two samples are drawn from subpopulations with equal variances.
      static boolean homoscedasticTTest​(double[] sample1, double[] sample2, double alpha)
      Performs a two-sided t-test evaluating the null hypothesis that sample1 and sample2 are drawn from populations with the same mean, with significance level alpha, assuming that the subpopulation variances are equal.
      static double homoscedasticTTest​(StatisticalSummary sampleStats1, StatisticalSummary sampleStats2)
      Returns the observed significance level, or p-value, associated with a two-sample, two-tailed t-test comparing the means of the datasets described by two StatisticalSummary instances, under the hypothesis of equal subpopulation variances.
      static double kolmogorovSmirnovStatistic​(double[] x, double[] y)
      Computes the two-sample Kolmogorov-Smirnov test statistic, \(D_{n,m}=\sup_x |F_n(x)-F_m(x)|\) where \(n\) is the length of x, \(m\) is the length of y, \(F_n\) is the empirical distribution that puts mass \(1/n\) at each of the values in x and \(F_m\) is the empirical distribution of the y values.
      static double kolmogorovSmirnovStatistic​(RealDistribution dist, double[] data)
      Computes the one-sample Kolmogorov-Smirnov test statistic, \(D_n=\sup_x |F_n(x)-F(x)|\) where \(F\) is the distribution (cdf) function associated with distribution, \(n\) is the length of data and \(F_n\) is the empirical distribution that puts mass \(1/n\) at each of the values in data.
      static double kolmogorovSmirnovTest​(double[] x, double[] y)
      Computes the p-value, or observed significance level, of a two-sample Kolmogorov-Smirnov test evaluating the null hypothesis that x and y are samples drawn from the same probability distribution.
      static double kolmogorovSmirnovTest​(double[] x, double[] y, boolean strict)
      Computes the p-value, or observed significance level, of a two-sample Kolmogorov-Smirnov test evaluating the null hypothesis that x and y are samples drawn from the same probability distribution.
      static double kolmogorovSmirnovTest​(RealDistribution dist, double[] data)
      Computes the p-value, or observed significance level, of a one-sample Kolmogorov-Smirnov test evaluating the null hypothesis that data conforms to distribution.
      static double kolmogorovSmirnovTest​(RealDistribution dist, double[] data, boolean strict)
      Computes the p-value, or observed significance level, of a one-sample Kolmogorov-Smirnov test evaluating the null hypothesis that data conforms to distribution.
      static boolean kolmogorovSmirnovTest​(RealDistribution dist, double[] data, double alpha)
      Performs a Kolmogorov-Smirnov test evaluating the null hypothesis that data conforms to distribution.
      static double oneWayAnovaFValue​(Collection<double[]> categoryData)
      Computes the ANOVA F-value for a collection of double[] arrays.
      static double oneWayAnovaPValue​(Collection<double[]> categoryData)
      Computes the ANOVA P-value for a collection of double[] arrays.
      static boolean oneWayAnovaTest​(Collection<double[]> categoryData, double alpha)
      Performs an ANOVA test, evaluating the null hypothesis that there is no difference among the means of the data categories.
      static double pairedT​(double[] sample1, double[] sample2)
      Computes a paired, 2-sample t-statistic based on the data in the input arrays.
      static double pairedTTest​(double[] sample1, double[] sample2)
      Returns the observed significance level, or p-value, associated with a paired, two-sample, two-tailed t-test based on the data in the input arrays.
      static boolean pairedTTest​(double[] sample1, double[] sample2, double alpha)
      Performs a paired t-test evaluating the null hypothesis that the mean of the paired differences between sample1 and sample2 is 0 in favor of the two-sided alternative that the mean paired difference is not equal to 0, with significance level alpha.
      static double rootLogLikelihoodRatio​(long k11, long k12, long k21, long k22)
      Calculates the root log-likelihood ratio for 2 state Datasets.
      static double t​(double[] sample1, double[] sample2)
      Computes a 2-sample t statistic, without the hypothesis of equal subpopulation variances.
      static double t​(double mu, double[] observed)
      Computes a t statistic given observed values and a comparison constant.
      static double t​(double mu, StatisticalSummary sampleStats)
      Computes a t statistic to use in comparing the mean of the dataset described by sampleStats to mu.
      static double t​(StatisticalSummary sampleStats1, StatisticalSummary sampleStats2)
      Computes a 2-sample t statistic, comparing the means of the datasets described by two StatisticalSummary instances, without the assumption of equal subpopulation variances.
      static double tTest​(double[] sample1, double[] sample2)
      Returns the observed significance level, or p-value, associated with a two-sample, two-tailed t-test comparing the means of the input arrays.
      static boolean tTest​(double[] sample1, double[] sample2, double alpha)
      Performs a two-sided t-test evaluating the null hypothesis that sample1 and sample2 are drawn from populations with the same mean, with significance level alpha.
      static double tTest​(double mu, double[] sample)
      Returns the observed significance level, or p-value, associated with a one-sample, two-tailed t-test comparing the mean of the input array with the constant mu.
      static boolean tTest​(double mu, double[] sample, double alpha)
      Performs a two-sided t-test evaluating the null hypothesis that the mean of the population from which sample is drawn equals mu.
      static double tTest​(double mu, StatisticalSummary sampleStats)
      Returns the observed significance level, or p-value, associated with a one-sample, two-tailed t-test comparing the mean of the dataset described by sampleStats with the constant mu.
      static boolean tTest​(double mu, StatisticalSummary sampleStats, double alpha)
      Performs a two-sided t-test evaluating the null hypothesis that the mean of the population from which the dataset described by stats is drawn equals mu.
      static double tTest​(StatisticalSummary sampleStats1, StatisticalSummary sampleStats2)
      Returns the observed significance level, or p-value, associated with a two-sample, two-tailed t-test comparing the means of the datasets described by two StatisticalSummary instances.
      static boolean tTest​(StatisticalSummary sampleStats1, StatisticalSummary sampleStats2, double alpha)
      Performs a two-sided t-test evaluating the null hypothesis that sampleStats1 and sampleStats2 describe datasets drawn from populations with the same mean, with significance level alpha.
    • Method Detail

      • homoscedasticT

        public static double homoscedasticT​(double[] sample1,
                                            double[] sample2)
                                     throws MathIllegalArgumentException,
                                            NullArgumentException
        Computes a 2-sample t statistic, under the hypothesis of equal subpopulation variances. To compute a t-statistic without the equal variances hypothesis, use t(double[], double[]).

        This statistic can be used to perform a (homoscedastic) two-sample t-test to compare sample means.

        The t-statistic is

           t = (m1 - m2) / (sqrt(1/n1 +1/n2) sqrt(var))

        where n1 is the size of first sample; n2 is the size of second sample; m1 is the mean of first sample; m2 is the mean of second sample and var is the pooled variance estimate:

        var = sqrt(((n1 - 1)var1 + (n2 - 1)var2) / ((n1-1) + (n2-1)))

        with var1 the variance of the first sample and var2 the variance of the second sample.

        Preconditions:

        • The observed array lengths must both be at least 2.
        Parameters:
        sample1 - array of sample data values
        sample2 - array of sample data values
        Returns:
        t statistic
        Throws:
        NullArgumentException - if the arrays are null
        MathIllegalArgumentException - if the length of the arrays is < 2
      • homoscedasticT

        public static double homoscedasticT​(StatisticalSummary sampleStats1,
                                            StatisticalSummary sampleStats2)
                                     throws MathIllegalArgumentException,
                                            NullArgumentException
        Computes a 2-sample t statistic, comparing the means of the datasets described by two StatisticalSummary instances, under the assumption of equal subpopulation variances. To compute a t-statistic without the equal variances assumption, use t(StatisticalSummary, StatisticalSummary).

        This statistic can be used to perform a (homoscedastic) two-sample t-test to compare sample means.

        The t-statistic returned is

           t = (m1 - m2) / (sqrt(1/n1 +1/n2) sqrt(var))

        where n1 is the size of first sample; n2 is the size of second sample; m1 is the mean of first sample; m2 is the mean of second sample and var is the pooled variance estimate:

        var = sqrt(((n1 - 1)var1 + (n2 - 1)var2) / ((n1-1) + (n2-1)))

        with var1 the variance of the first sample and var2 the variance of the second sample.

        Preconditions:

        • The datasets described by the two Univariates must each contain at least 2 observations.
        Parameters:
        sampleStats1 - StatisticalSummary describing data from the first sample
        sampleStats2 - StatisticalSummary describing data from the second sample
        Returns:
        t statistic
        Throws:
        NullArgumentException - if the sample statistics are null
        MathIllegalArgumentException - if the number of samples is < 2
      • homoscedasticTTest

        public static boolean homoscedasticTTest​(double[] sample1,
                                                 double[] sample2,
                                                 double alpha)
                                          throws MathIllegalArgumentException,
                                                 NullArgumentException,
                                                 MathIllegalStateException
        Performs a two-sided t-test evaluating the null hypothesis that sample1 and sample2 are drawn from populations with the same mean, with significance level alpha, assuming that the subpopulation variances are equal. Use tTest(double[], double[], double) to perform the test without the assumption of equal variances.

        Returns true iff the null hypothesis that the means are equal can be rejected with confidence 1 - alpha. To perform a 1-sided test, use alpha * 2. To perform the test without the assumption of equal subpopulation variances, use tTest(double[], double[], double).

        A pooled variance estimate is used to compute the t-statistic. See t(double[], double[]) for the formula. The sum of the sample sizes minus 2 is used as the degrees of freedom.

        Examples:

        1. To test the (2-sided) hypothesis mean 1 = mean 2 at the 95% level, use
          tTest(sample1, sample2, 0.05).
        2. To test the (one-sided) hypothesis mean 1 < mean 2, at the 99% level, first verify that the measured mean of sample 1 is less than the mean of sample 2 and then use
          tTest(sample1, sample2, 0.02)

        Usage Note:
        The validity of the test depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The observed array lengths must both be at least 2.
        • 0 < alpha < 0.5
        Parameters:
        sample1 - array of sample data values
        sample2 - array of sample data values
        alpha - significance level of the test
        Returns:
        true if the null hypothesis can be rejected with confidence 1 - alpha
        Throws:
        NullArgumentException - if the arrays are null
        MathIllegalArgumentException - if the length of the arrays is < 2
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if an error occurs computing the p-value
      • homoscedasticTTest

        public static double homoscedasticTTest​(double[] sample1,
                                                double[] sample2)
                                         throws MathIllegalArgumentException,
                                                NullArgumentException,
                                                MathIllegalStateException
        Returns the observed significance level, or p-value, associated with a two-sample, two-tailed t-test comparing the means of the input arrays, under the assumption that the two samples are drawn from subpopulations with equal variances. To perform the test without the equal variances assumption, use tTest(double[], double[]).

        The number returned is the smallest significance level at which one can reject the null hypothesis that the two means are equal in favor of the two-sided alternative that they are different. For a one-sided test, divide the returned value by 2.

        A pooled variance estimate is used to compute the t-statistic. See homoscedasticT(double[], double[]). The sum of the sample sizes minus 2 is used as the degrees of freedom.

        Usage Note:
        The validity of the p-value depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The observed array lengths must both be at least 2.
        Parameters:
        sample1 - array of sample data values
        sample2 - array of sample data values
        Returns:
        p-value for t-test
        Throws:
        NullArgumentException - if the arrays are null
        MathIllegalArgumentException - if the length of the arrays is < 2
        MathIllegalStateException - if an error occurs computing the p-value
      • homoscedasticTTest

        public static double homoscedasticTTest​(StatisticalSummary sampleStats1,
                                                StatisticalSummary sampleStats2)
                                         throws MathIllegalArgumentException,
                                                NullArgumentException,
                                                MathIllegalStateException
        Returns the observed significance level, or p-value, associated with a two-sample, two-tailed t-test comparing the means of the datasets described by two StatisticalSummary instances, under the hypothesis of equal subpopulation variances. To perform a test without the equal variances assumption, use tTest(StatisticalSummary, StatisticalSummary).

        The number returned is the smallest significance level at which one can reject the null hypothesis that the two means are equal in favor of the two-sided alternative that they are different. For a one-sided test, divide the returned value by 2.

        See homoscedasticT(double[], double[]) for the formula used to compute the t-statistic. The sum of the sample sizes minus 2 is used as the degrees of freedom.

        Usage Note:
        The validity of the p-value depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The datasets described by the two Univariates must each contain at least 2 observations.
        Parameters:
        sampleStats1 - StatisticalSummary describing data from the first sample
        sampleStats2 - StatisticalSummary describing data from the second sample
        Returns:
        p-value for t-test
        Throws:
        NullArgumentException - if the sample statistics are null
        MathIllegalArgumentException - if the number of samples is < 2
        MathIllegalStateException - if an error occurs computing the p-value
      • pairedT

        public static double pairedT​(double[] sample1,
                                     double[] sample2)
                              throws MathIllegalArgumentException,
                                     NullArgumentException
        Computes a paired, 2-sample t-statistic based on the data in the input arrays. The t-statistic returned is equivalent to what would be returned by computing the one-sample t-statistic t(double, double[]), with mu = 0 and the sample array consisting of the (signed) differences between corresponding entries in sample1 and sample2.

        Preconditions:

        • The input arrays must have the same length and their common length must be at least 2.
        Parameters:
        sample1 - array of sample data values
        sample2 - array of sample data values
        Returns:
        t statistic
        Throws:
        NullArgumentException - if the arrays are null
        MathIllegalArgumentException - if the arrays are empty
        MathIllegalArgumentException - if the length of the arrays is not equal
        MathIllegalArgumentException - if the length of the arrays is < 2
      • pairedTTest

        public static boolean pairedTTest​(double[] sample1,
                                          double[] sample2,
                                          double alpha)
                                   throws MathIllegalArgumentException,
                                          NullArgumentException,
                                          MathIllegalStateException
        Performs a paired t-test evaluating the null hypothesis that the mean of the paired differences between sample1 and sample2 is 0 in favor of the two-sided alternative that the mean paired difference is not equal to 0, with significance level alpha.

        Returns true iff the null hypothesis can be rejected with confidence 1 - alpha. To perform a 1-sided test, use alpha * 2

        Usage Note:
        The validity of the test depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The input array lengths must be the same and their common length must be at least 2.
        • 0 < alpha < 0.5
        Parameters:
        sample1 - array of sample data values
        sample2 - array of sample data values
        alpha - significance level of the test
        Returns:
        true if the null hypothesis can be rejected with confidence 1 - alpha
        Throws:
        NullArgumentException - if the arrays are null
        MathIllegalArgumentException - if the arrays are empty
        MathIllegalArgumentException - if the length of the arrays is not equal
        MathIllegalArgumentException - if the length of the arrays is < 2
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if an error occurs computing the p-value
      • pairedTTest

        public static double pairedTTest​(double[] sample1,
                                         double[] sample2)
                                  throws MathIllegalArgumentException,
                                         NullArgumentException,
                                         MathIllegalStateException
        Returns the observed significance level, or p-value, associated with a paired, two-sample, two-tailed t-test based on the data in the input arrays.

        The number returned is the smallest significance level at which one can reject the null hypothesis that the mean of the paired differences is 0 in favor of the two-sided alternative that the mean paired difference is not equal to 0. For a one-sided test, divide the returned value by 2.

        This test is equivalent to a one-sample t-test computed using tTest(double, double[]) with mu = 0 and the sample array consisting of the signed differences between corresponding elements of sample1 and sample2.

        Usage Note:
        The validity of the p-value depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The input array lengths must be the same and their common length must be at least 2.
        Parameters:
        sample1 - array of sample data values
        sample2 - array of sample data values
        Returns:
        p-value for t-test
        Throws:
        NullArgumentException - if the arrays are null
        MathIllegalArgumentException - if the arrays are empty
        MathIllegalArgumentException - if the length of the arrays is not equal
        MathIllegalArgumentException - if the length of the arrays is < 2
        MathIllegalStateException - if an error occurs computing the p-value
      • t

        public static double t​(double[] sample1,
                               double[] sample2)
                        throws MathIllegalArgumentException,
                               NullArgumentException
        Computes a 2-sample t statistic, without the hypothesis of equal subpopulation variances. To compute a t-statistic assuming equal variances, use homoscedasticT(double[], double[]).

        This statistic can be used to perform a two-sample t-test to compare sample means.

        The t-statistic is

           t = (m1 - m2) / sqrt(var1/n1 + var2/n2)

        where n1 is the size of the first sample n2 is the size of the second sample; m1 is the mean of the first sample; m2 is the mean of the second sample; var1 is the variance of the first sample; var2 is the variance of the second sample;

        Preconditions:

        • The observed array lengths must both be at least 2.
        Parameters:
        sample1 - array of sample data values
        sample2 - array of sample data values
        Returns:
        t statistic
        Throws:
        NullArgumentException - if the arrays are null
        MathIllegalArgumentException - if the length of the arrays is < 2
      • t

        public static double t​(StatisticalSummary sampleStats1,
                               StatisticalSummary sampleStats2)
                        throws MathIllegalArgumentException,
                               NullArgumentException
        Computes a 2-sample t statistic, comparing the means of the datasets described by two StatisticalSummary instances, without the assumption of equal subpopulation variances. Use homoscedasticT(StatisticalSummary, StatisticalSummary) to compute a t-statistic under the equal variances assumption.

        This statistic can be used to perform a two-sample t-test to compare sample means.

        The returned t-statistic is

           t = (m1 - m2) / sqrt(var1/n1 + var2/n2)

        where n1 is the size of the first sample; n2 is the size of the second sample; m1 is the mean of the first sample; m2 is the mean of the second sample var1 is the variance of the first sample; var2 is the variance of the second sample

        Preconditions:

        • The datasets described by the two Univariates must each contain at least 2 observations.
        Parameters:
        sampleStats1 - StatisticalSummary describing data from the first sample
        sampleStats2 - StatisticalSummary describing data from the second sample
        Returns:
        t statistic
        Throws:
        NullArgumentException - if the sample statistics are null
        MathIllegalArgumentException - if the number of samples is < 2
      • tTest

        public static boolean tTest​(double mu,
                                    double[] sample,
                                    double alpha)
                             throws MathIllegalArgumentException,
                                    NullArgumentException,
                                    MathIllegalStateException
        Performs a two-sided t-test evaluating the null hypothesis that the mean of the population from which sample is drawn equals mu.

        Returns true iff the null hypothesis can be rejected with confidence 1 - alpha. To perform a 1-sided test, use alpha * 2

        Examples:

        1. To test the (2-sided) hypothesis sample mean = mu at the 95% level, use
          tTest(mu, sample, 0.05)
        2. To test the (one-sided) hypothesis sample mean < mu at the 99% level, first verify that the measured sample mean is less than mu and then use
          tTest(mu, sample, 0.02)

        Usage Note:
        The validity of the test depends on the assumptions of the one-sample parametric t-test procedure, as discussed here

        Preconditions:

        • The observed array length must be at least 2.
        Parameters:
        mu - constant value to compare sample mean against
        sample - array of sample data values
        alpha - significance level of the test
        Returns:
        p-value
        Throws:
        NullArgumentException - if the sample array is null
        MathIllegalArgumentException - if the length of the array is < 2
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if an error computing the p-value
      • tTest

        public static double tTest​(double mu,
                                   double[] sample)
                            throws MathIllegalArgumentException,
                                   NullArgumentException,
                                   MathIllegalStateException
        Returns the observed significance level, or p-value, associated with a one-sample, two-tailed t-test comparing the mean of the input array with the constant mu.

        The number returned is the smallest significance level at which one can reject the null hypothesis that the mean equals mu in favor of the two-sided alternative that the mean is different from mu. For a one-sided test, divide the returned value by 2.

        Usage Note:
        The validity of the test depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The observed array length must be at least 2.
        Parameters:
        mu - constant value to compare sample mean against
        sample - array of sample data values
        Returns:
        p-value
        Throws:
        NullArgumentException - if the sample array is null
        MathIllegalArgumentException - if the length of the array is < 2
        MathIllegalStateException - if an error occurs computing the p-value
      • tTest

        public static boolean tTest​(double mu,
                                    StatisticalSummary sampleStats,
                                    double alpha)
                             throws MathIllegalArgumentException,
                                    NullArgumentException,
                                    MathIllegalStateException
        Performs a two-sided t-test evaluating the null hypothesis that the mean of the population from which the dataset described by stats is drawn equals mu.

        Returns true iff the null hypothesis can be rejected with confidence 1 - alpha. To perform a 1-sided test, use alpha * 2.

        Examples:

        1. To test the (2-sided) hypothesis sample mean = mu at the 95% level, use
          tTest(mu, sampleStats, 0.05)
        2. To test the (one-sided) hypothesis sample mean < mu at the 99% level, first verify that the measured sample mean is less than mu and then use
          tTest(mu, sampleStats, 0.02)

        Usage Note:
        The validity of the test depends on the assumptions of the one-sample parametric t-test procedure, as discussed here

        Preconditions:

        • The sample must include at least 2 observations.
        Parameters:
        mu - constant value to compare sample mean against
        sampleStats - StatisticalSummary describing sample data values
        alpha - significance level of the test
        Returns:
        p-value
        Throws:
        NullArgumentException - if sampleStats is null
        MathIllegalArgumentException - if the number of samples is < 2
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if an error occurs computing the p-value
      • tTest

        public static double tTest​(double mu,
                                   StatisticalSummary sampleStats)
                            throws MathIllegalArgumentException,
                                   NullArgumentException,
                                   MathIllegalStateException
        Returns the observed significance level, or p-value, associated with a one-sample, two-tailed t-test comparing the mean of the dataset described by sampleStats with the constant mu.

        The number returned is the smallest significance level at which one can reject the null hypothesis that the mean equals mu in favor of the two-sided alternative that the mean is different from mu. For a one-sided test, divide the returned value by 2.

        Usage Note:
        The validity of the test depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The sample must contain at least 2 observations.
        Parameters:
        mu - constant value to compare sample mean against
        sampleStats - StatisticalSummary describing sample data
        Returns:
        p-value
        Throws:
        NullArgumentException - if sampleStats is null
        MathIllegalArgumentException - if the number of samples is < 2
        MathIllegalStateException - if an error occurs computing the p-value
      • tTest

        public static boolean tTest​(double[] sample1,
                                    double[] sample2,
                                    double alpha)
                             throws MathIllegalArgumentException,
                                    NullArgumentException,
                                    MathIllegalStateException
        Performs a two-sided t-test evaluating the null hypothesis that sample1 and sample2 are drawn from populations with the same mean, with significance level alpha. This test does not assume that the subpopulation variances are equal. To perform the test assuming equal variances, use homoscedasticTTest(double[], double[], double).

        Returns true iff the null hypothesis that the means are equal can be rejected with confidence 1 - alpha. To perform a 1-sided test, use alpha * 2

        See t(double[], double[]) for the formula used to compute the t-statistic. Degrees of freedom are approximated using the Welch-Satterthwaite approximation.

        Examples:

        1. To test the (2-sided) hypothesis mean 1 = mean 2 at the 95% level, use
          tTest(sample1, sample2, 0.05).
        2. To test the (one-sided) hypothesis mean 1 < mean 2 , at the 99% level, first verify that the measured mean of sample 1 is less than the mean of sample 2 and then use
          tTest(sample1, sample2, 0.02)

        Usage Note:
        The validity of the test depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The observed array lengths must both be at least 2.
        • 0 < alpha < 0.5
        Parameters:
        sample1 - array of sample data values
        sample2 - array of sample data values
        alpha - significance level of the test
        Returns:
        true if the null hypothesis can be rejected with confidence 1 - alpha
        Throws:
        NullArgumentException - if the arrays are null
        MathIllegalArgumentException - if the length of the arrays is < 2
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if an error occurs computing the p-value
      • tTest

        public static double tTest​(double[] sample1,
                                   double[] sample2)
                            throws MathIllegalArgumentException,
                                   NullArgumentException,
                                   MathIllegalStateException
        Returns the observed significance level, or p-value, associated with a two-sample, two-tailed t-test comparing the means of the input arrays.

        The number returned is the smallest significance level at which one can reject the null hypothesis that the two means are equal in favor of the two-sided alternative that they are different. For a one-sided test, divide the returned value by 2.

        The test does not assume that the underlying popuation variances are equal and it uses approximated degrees of freedom computed from the sample data to compute the p-value. The t-statistic used is as defined in t(double[], double[]) and the Welch-Satterthwaite approximation to the degrees of freedom is used, as described here. To perform the test under the assumption of equal subpopulation variances, use homoscedasticTTest(double[], double[]).

        Usage Note:
        The validity of the p-value depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The observed array lengths must both be at least 2.
        Parameters:
        sample1 - array of sample data values
        sample2 - array of sample data values
        Returns:
        p-value for t-test
        Throws:
        NullArgumentException - if the arrays are null
        MathIllegalArgumentException - if the length of the arrays is < 2
        MathIllegalStateException - if an error occurs computing the p-value
      • tTest

        public static boolean tTest​(StatisticalSummary sampleStats1,
                                    StatisticalSummary sampleStats2,
                                    double alpha)
                             throws MathIllegalArgumentException,
                                    NullArgumentException,
                                    MathIllegalStateException
        Performs a two-sided t-test evaluating the null hypothesis that sampleStats1 and sampleStats2 describe datasets drawn from populations with the same mean, with significance level alpha. This test does not assume that the subpopulation variances are equal. To perform the test under the equal variances assumption, use homoscedasticTTest(StatisticalSummary, StatisticalSummary).

        Returns true iff the null hypothesis that the means are equal can be rejected with confidence 1 - alpha. To perform a 1-sided test, use alpha * 2

        See t(double[], double[]) for the formula used to compute the t-statistic. Degrees of freedom are approximated using the Welch-Satterthwaite approximation.

        Examples:

        1. To test the (2-sided) hypothesis mean 1 = mean 2 at the 95%, use
          tTest(sampleStats1, sampleStats2, 0.05)
        2. To test the (one-sided) hypothesis mean 1 < mean 2 at the 99% level, first verify that the measured mean of sample 1 is less than the mean of sample 2 and then use
          tTest(sampleStats1, sampleStats2, 0.02)

        Usage Note:
        The validity of the test depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The datasets described by the two Univariates must each contain at least 2 observations.
        • 0 < alpha < 0.5
        Parameters:
        sampleStats1 - StatisticalSummary describing sample data values
        sampleStats2 - StatisticalSummary describing sample data values
        alpha - significance level of the test
        Returns:
        true if the null hypothesis can be rejected with confidence 1 - alpha
        Throws:
        NullArgumentException - if the sample statistics are null
        MathIllegalArgumentException - if the number of samples is < 2
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if an error occurs computing the p-value
      • tTest

        public static double tTest​(StatisticalSummary sampleStats1,
                                   StatisticalSummary sampleStats2)
                            throws MathIllegalArgumentException,
                                   NullArgumentException,
                                   MathIllegalStateException
        Returns the observed significance level, or p-value, associated with a two-sample, two-tailed t-test comparing the means of the datasets described by two StatisticalSummary instances.

        The number returned is the smallest significance level at which one can reject the null hypothesis that the two means are equal in favor of the two-sided alternative that they are different. For a one-sided test, divide the returned value by 2.

        The test does not assume that the underlying population variances are equal and it uses approximated degrees of freedom computed from the sample data to compute the p-value. To perform the test assuming equal variances, use homoscedasticTTest(StatisticalSummary, StatisticalSummary).

        Usage Note:
        The validity of the p-value depends on the assumptions of the parametric t-test procedure, as discussed here

        Preconditions:

        • The datasets described by the two Univariates must each contain at least 2 observations.
        Parameters:
        sampleStats1 - StatisticalSummary describing data from the first sample
        sampleStats2 - StatisticalSummary describing data from the second sample
        Returns:
        p-value for t-test
        Throws:
        NullArgumentException - if the sample statistics are null
        MathIllegalArgumentException - if the number of samples is < 2
        MathIllegalStateException - if an error occurs computing the p-value
      • chiSquare

        public static double chiSquare​(double[] expected,
                                       long[] observed)
                                throws MathIllegalArgumentException
        Computes the Chi-Square statistic comparing observed and expected frequency counts.

        This statistic can be used to perform a Chi-Square test evaluating the null hypothesis that the observed counts follow the expected distribution.

        Preconditions:

        • Expected counts must all be positive.
        • Observed counts must all be ≥ 0.
        • The observed and expected arrays must have the same length and their common length must be at least 2.

        If any of the preconditions are not met, an IllegalArgumentException is thrown.

        Note: This implementation rescales the expected array if necessary to ensure that the sum of the expected and observed counts are equal.

        Parameters:
        observed - array of observed frequency counts
        expected - array of expected frequency counts
        Returns:
        chiSquare test statistic
        Throws:
        MathIllegalArgumentException - if observed has negative entries
        MathIllegalArgumentException - if expected has entries that are not strictly positive
        MathIllegalArgumentException - if the arrays length is less than 2
      • chiSquare

        public static double chiSquare​(long[][] counts)
                                throws MathIllegalArgumentException,
                                       NullArgumentException
        Computes the Chi-Square statistic associated with a chi-square test of independence based on the input counts array, viewed as a two-way table.

        The rows of the 2-way table are count[0], ... , count[count.length - 1]

        Preconditions:

        • All counts must be ≥ 0.
        • The count array must be rectangular (i.e. all count[i] subarrays must have the same length).
        • The 2-way table represented by counts must have at least 2 columns and at least 2 rows.

        If any of the preconditions are not met, an IllegalArgumentException is thrown.

        Parameters:
        counts - array representation of 2-way table
        Returns:
        chiSquare test statistic
        Throws:
        NullArgumentException - if the array is null
        MathIllegalArgumentException - if the array is not rectangular
        MathIllegalArgumentException - if counts has negative entries
      • chiSquareTest

        public static boolean chiSquareTest​(double[] expected,
                                            long[] observed,
                                            double alpha)
                                     throws MathIllegalArgumentException,
                                            MathIllegalStateException
        Performs a Chi-square goodness of fit test evaluating the null hypothesis that the observed counts conform to the frequency distribution described by the expected counts, with significance level alpha. Returns true iff the null hypothesis can be rejected with 100 * (1 - alpha) percent confidence.

        Example:
        To test the hypothesis that observed follows expected at the 99% level, use chiSquareTest(expected, observed, 0.01)

        Preconditions:

        • Expected counts must all be positive.
        • Observed counts must all be ≥ 0.
        • The observed and expected arrays must have the same length and their common length must be at least 2.
        • 0 < alpha < 0.5

        If any of the preconditions are not met, an IllegalArgumentException is thrown.

        Note: This implementation rescales the expected array if necessary to ensure that the sum of the expected and observed counts are equal.

        Parameters:
        observed - array of observed frequency counts
        expected - array of expected frequency counts
        alpha - significance level of the test
        Returns:
        true iff null hypothesis can be rejected with confidence 1 - alpha
        Throws:
        MathIllegalArgumentException - if observed has negative entries
        MathIllegalArgumentException - if expected has entries that are not strictly positive
        MathIllegalArgumentException - if the arrays length is less than 2
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if an error occurs computing the p-value
      • chiSquareTest

        public static double chiSquareTest​(double[] expected,
                                           long[] observed)
                                    throws MathIllegalArgumentException,
                                           MathIllegalStateException
        Returns the observed significance level, or p-value, associated with a Chi-square goodness of fit test comparing the observed frequency counts to those in the expected array.

        The number returned is the smallest significance level at which one can reject the null hypothesis that the observed counts conform to the frequency distribution described by the expected counts.

        Preconditions:

        • Expected counts must all be positive.
        • Observed counts must all be ≥ 0.
        • The observed and expected arrays must have the same length and their common length must be at least 2.

        If any of the preconditions are not met, an IllegalArgumentException is thrown.

        Note: This implementation rescales the expected array if necessary to ensure that the sum of the expected and observed counts are equal.

        Parameters:
        observed - array of observed frequency counts
        expected - array of expected frequency counts
        Returns:
        p-value
        Throws:
        MathIllegalArgumentException - if observed has negative entries
        MathIllegalArgumentException - if expected has entries that are not strictly positive
        MathIllegalArgumentException - if the arrays length is less than 2
        MathIllegalStateException - if an error occurs computing the p-value
      • chiSquareTest

        public static boolean chiSquareTest​(long[][] counts,
                                            double alpha)
                                     throws MathIllegalArgumentException,
                                            NullArgumentException,
                                            MathIllegalStateException
        Performs a chi-square test of independence evaluating the null hypothesis that the classifications represented by the counts in the columns of the input 2-way table are independent of the rows, with significance level alpha. Returns true iff the null hypothesis can be rejected with 100 * (1 - alpha) percent confidence.

        The rows of the 2-way table are count[0], ... , count[count.length - 1]

        Example:
        To test the null hypothesis that the counts in count[0], ... , count[count.length - 1] all correspond to the same underlying probability distribution at the 99% level, use chiSquareTest(counts, 0.01).

        Preconditions:

        • All counts must be ≥ 0.
        • The count array must be rectangular (i.e. all count[i] subarrays must have the same length).
        • The 2-way table represented by counts must have at least 2 columns and at least 2 rows.

        If any of the preconditions are not met, an IllegalArgumentException is thrown.

        Parameters:
        counts - array representation of 2-way table
        alpha - significance level of the test
        Returns:
        true iff null hypothesis can be rejected with confidence 1 - alpha
        Throws:
        NullArgumentException - if the array is null
        MathIllegalArgumentException - if the array is not rectangular
        MathIllegalArgumentException - if counts has any negative entries
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if an error occurs computing the p-value
      • chiSquareDataSetsComparison

        public static double chiSquareDataSetsComparison​(long[] observed1,
                                                         long[] observed2)
                                                  throws MathIllegalArgumentException
        Computes a Chi-Square two sample test statistic comparing bin frequency counts in observed1 and observed2.

        The sums of frequency counts in the two samples are not required to be the same. The formula used to compute the test statistic is

        ∑[(K * observed1[i] - observed2[i]/K)2 / (observed1[i] + observed2[i])]

        where

        K = √[∑(observed2 / ∑(observed1)]

        This statistic can be used to perform a Chi-Square test evaluating the null hypothesis that both observed counts follow the same distribution.

        Preconditions:

        • Observed counts must be non-negative.
        • Observed counts for a specific bin must not both be zero.
        • Observed counts for a specific sample must not all be 0.
        • The arrays observed1 and observed2 must have the same length and their common length must be at least 2.

        If any of the preconditions are not met, an IllegalArgumentException is thrown.

        Parameters:
        observed1 - array of observed frequency counts of the first data set
        observed2 - array of observed frequency counts of the second data set
        Returns:
        chiSquare test statistic
        Throws:
        MathIllegalArgumentException - the the length of the arrays does not match
        MathIllegalArgumentException - if any entries in observed1 or observed2 are negative
        MathIllegalArgumentException - if either all counts of observed1 or observed2 are zero, or if the count at some index is zero for both arrays
      • chiSquareTestDataSetsComparison

        public static double chiSquareTestDataSetsComparison​(long[] observed1,
                                                             long[] observed2)
                                                      throws MathIllegalArgumentException,
                                                             MathIllegalStateException
        Returns the observed significance level, or p-value, associated with a Chi-Square two sample test comparing bin frequency counts in observed1 and observed2.

        The number returned is the smallest significance level at which one can reject the null hypothesis that the observed counts conform to the same distribution.

        See chiSquareDataSetsComparison(long[], long[]) for details on the formula used to compute the test statistic. The degrees of of freedom used to perform the test is one less than the common length of the input observed count arrays.

        Preconditions:

        • Observed counts must be non-negative.
        • Observed counts for a specific bin must not both be zero.
        • Observed counts for a specific sample must not all be 0.
        • The arrays observed1 and observed2 must have the same length and their common length must be at least 2.

        If any of the preconditions are not met, an IllegalArgumentException is thrown.

        Parameters:
        observed1 - array of observed frequency counts of the first data set
        observed2 - array of observed frequency counts of the second data set
        Returns:
        p-value
        Throws:
        MathIllegalArgumentException - the the length of the arrays does not match
        MathIllegalArgumentException - if any entries in observed1 or observed2 are negative
        MathIllegalArgumentException - if either all counts of observed1 or observed2 are zero, or if the count at the same index is zero for both arrays
        MathIllegalStateException - if an error occurs computing the p-value
      • chiSquareTestDataSetsComparison

        public static boolean chiSquareTestDataSetsComparison​(long[] observed1,
                                                              long[] observed2,
                                                              double alpha)
                                                       throws MathIllegalArgumentException,
                                                              MathIllegalStateException
        Performs a Chi-Square two sample test comparing two binned data sets. The test evaluates the null hypothesis that the two lists of observed counts conform to the same frequency distribution, with significance level alpha. Returns true iff the null hypothesis can be rejected with 100 * (1 - alpha) percent confidence.

        See chiSquareDataSetsComparison(long[], long[]) for details on the formula used to compute the Chisquare statistic used in the test. The degrees of of freedom used to perform the test is one less than the common length of the input observed count arrays.

        Preconditions:

        • Observed counts must be non-negative.
        • Observed counts for a specific bin must not both be zero.
        • Observed counts for a specific sample must not all be 0.
        • The arrays observed1 and observed2 must have the same length and their common length must be at least 2.
        • 0 < alpha < 0.5

        If any of the preconditions are not met, an IllegalArgumentException is thrown.

        Parameters:
        observed1 - array of observed frequency counts of the first data set
        observed2 - array of observed frequency counts of the second data set
        alpha - significance level of the test
        Returns:
        true iff null hypothesis can be rejected with confidence 1 - alpha
        Throws:
        MathIllegalArgumentException - the the length of the arrays does not match
        MathIllegalArgumentException - if any entries in observed1 or observed2 are negative
        MathIllegalArgumentException - if either all counts of observed1 or observed2 are zero, or if the count at the same index is zero for both arrays
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if an error occurs performing the test
      • oneWayAnovaFValue

        public static double oneWayAnovaFValue​(Collection<double[]> categoryData)
                                        throws MathIllegalArgumentException,
                                               NullArgumentException
        Computes the ANOVA F-value for a collection of double[] arrays.

        Preconditions:

        • The categoryData Collection must contain double[] arrays.
        • There must be at least two double[] arrays in the categoryData collection and each of these arrays must contain at least two values.

        This implementation computes the F statistic using the definitional formula

           F = msbg/mswg

        where

          msbg = between group mean square
          mswg = within group mean square

        are as defined here

        Parameters:
        categoryData - Collection of double[] arrays each containing data for one category
        Returns:
        Fvalue
        Throws:
        NullArgumentException - if categoryData is null
        MathIllegalArgumentException - if the length of the categoryData array is less than 2 or a contained double[] array does not have at least two values
      • oneWayAnovaPValue

        public static double oneWayAnovaPValue​(Collection<double[]> categoryData)
                                        throws MathIllegalArgumentException,
                                               NullArgumentException,
                                               MathIllegalStateException
        Computes the ANOVA P-value for a collection of double[] arrays.

        Preconditions:

        • The categoryData Collection must contain double[] arrays.
        • There must be at least two double[] arrays in the categoryData collection and each of these arrays must contain at least two values.

        This implementation uses the Hipparchus F Distribution implementation to estimate the exact p-value, using the formula

           p = 1 - cumulativeProbability(F)

        where F is the F value and cumulativeProbability is the Hipparchus implementation of the F distribution.

        Parameters:
        categoryData - Collection of double[] arrays each containing data for one category
        Returns:
        Pvalue
        Throws:
        NullArgumentException - if categoryData is null
        MathIllegalArgumentException - if the length of the categoryData array is less than 2 or a contained double[] array does not have at least two values
        MathIllegalStateException - if the p-value can not be computed due to a convergence error
        MathIllegalStateException - if the maximum number of iterations is exceeded
      • oneWayAnovaTest

        public static boolean oneWayAnovaTest​(Collection<double[]> categoryData,
                                              double alpha)
                                       throws MathIllegalArgumentException,
                                              NullArgumentException,
                                              MathIllegalStateException
        Performs an ANOVA test, evaluating the null hypothesis that there is no difference among the means of the data categories.

        Preconditions:

        • The categoryData Collection must contain double[] arrays.
        • There must be at least two double[] arrays in the categoryData collection and each of these arrays must contain at least two values.
        • alpha must be strictly greater than 0 and less than or equal to 0.5.

        This implementation uses the Hipparchus F Distribution implementation to estimate the exact p-value, using the formula

           p = 1 - cumulativeProbability(F)

        where F is the F value and cumulativeProbability is the Hipparchus implementation of the F distribution.

        True is returned iff the estimated p-value is less than alpha.

        Parameters:
        categoryData - Collection of double[] arrays each containing data for one category
        alpha - significance level of the test
        Returns:
        true if the null hypothesis can be rejected with confidence 1 - alpha
        Throws:
        NullArgumentException - if categoryData is null
        MathIllegalArgumentException - if the length of the categoryData array is less than 2 or a contained double[] array does not have at least two values
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if the p-value can not be computed due to a convergence error
        MathIllegalStateException - if the maximum number of iterations is exceeded
      • g

        public static double g​(double[] expected,
                               long[] observed)
                        throws MathIllegalArgumentException
        Computes the G statistic for Goodness of Fit comparing observed and expected frequency counts.

        This statistic can be used to perform a G test (Log-Likelihood Ratio Test) evaluating the null hypothesis that the observed counts follow the expected distribution.

        Preconditions:

        • Expected counts must all be positive.
        • Observed counts must all be ≥ 0.
        • The observed and expected arrays must have the same length and their common length must be at least 2.

        If any of the preconditions are not met, a MathIllegalArgumentException is thrown.

        Note:This implementation rescales the expected array if necessary to ensure that the sum of the expected and observed counts are equal.

        Parameters:
        observed - array of observed frequency counts
        expected - array of expected frequency counts
        Returns:
        G-Test statistic
        Throws:
        MathIllegalArgumentException - if observed has negative entries
        MathIllegalArgumentException - if expected has entries that are not strictly positive
        MathIllegalArgumentException - if the array lengths do not match or are less than 2.
      • gTest

        public static double gTest​(double[] expected,
                                   long[] observed)
                            throws MathIllegalArgumentException,
                                   MathIllegalStateException
        Returns the observed significance level, or p-value, associated with a G-Test for goodness of fit comparing the observed frequency counts to those in the expected array.

        The number returned is the smallest significance level at which one can reject the null hypothesis that the observed counts conform to the frequency distribution described by the expected counts.

        The probability returned is the tail probability beyond g(expected, observed) in the ChiSquare distribution with degrees of freedom one less than the common length of expected and observed.

        Preconditions:

        • Expected counts must all be positive.
        • Observed counts must all be ≥ 0.
        • The observed and expected arrays must have the same length and their common length must be at least 2.

        If any of the preconditions are not met, a MathIllegalArgumentException is thrown.

        Note:This implementation rescales the expected array if necessary to ensure that the sum of the expected and observed counts are equal.

        Parameters:
        observed - array of observed frequency counts
        expected - array of expected frequency counts
        Returns:
        p-value
        Throws:
        MathIllegalArgumentException - if observed has negative entries
        MathIllegalArgumentException - if expected has entries that are not strictly positive
        MathIllegalArgumentException - if the array lengths do not match or are less than 2.
        MathIllegalStateException - if an error occurs computing the p-value.
      • gTestIntrinsic

        public static double gTestIntrinsic​(double[] expected,
                                            long[] observed)
                                     throws MathIllegalArgumentException,
                                            MathIllegalStateException
        Returns the intrinsic (Hardy-Weinberg proportions) p-Value, as described in p64-69 of McDonald, J.H. 2009. Handbook of Biological Statistics (2nd ed.). Sparky House Publishing, Baltimore, Maryland.

        The probability returned is the tail probability beyond g(expected, observed) in the ChiSquare distribution with degrees of freedom two less than the common length of expected and observed.

        Parameters:
        observed - array of observed frequency counts
        expected - array of expected frequency counts
        Returns:
        p-value
        Throws:
        MathIllegalArgumentException - if observed has negative entries
        MathIllegalArgumentException - expected has entries that are not strictly positive
        MathIllegalArgumentException - if the array lengths do not match or are less than 2.
        MathIllegalStateException - if an error occurs computing the p-value.
      • gTest

        public static boolean gTest​(double[] expected,
                                    long[] observed,
                                    double alpha)
                             throws MathIllegalArgumentException,
                                    MathIllegalStateException
        Performs a G-Test (Log-Likelihood Ratio Test) for goodness of fit evaluating the null hypothesis that the observed counts conform to the frequency distribution described by the expected counts, with significance level alpha. Returns true iff the null hypothesis can be rejected with 100 * (1 - alpha) percent confidence.

        Example:
        To test the hypothesis that observed follows expected at the 99% level, use

        gTest(expected, observed, 0.01)

        Returns true iff gTestGoodnessOfFitPValue(expected, observed) > alpha

        Preconditions:

        • Expected counts must all be positive.
        • Observed counts must all be ≥ 0.
        • The observed and expected arrays must have the same length and their common length must be at least 2.
        • 0 < alpha < 0.5

        If any of the preconditions are not met, a MathIllegalArgumentException is thrown.

        Note:This implementation rescales the expected array if necessary to ensure that the sum of the expected and observed counts are equal.

        Parameters:
        observed - array of observed frequency counts
        expected - array of expected frequency counts
        alpha - significance level of the test
        Returns:
        true iff null hypothesis can be rejected with confidence 1 - alpha
        Throws:
        MathIllegalArgumentException - if observed has negative entries
        MathIllegalArgumentException - if expected has entries that are not strictly positive
        MathIllegalArgumentException - if the array lengths do not match or are less than 2.
        MathIllegalStateException - if an error occurs computing the p-value.
        MathIllegalArgumentException - if alpha is not strictly greater than zero and less than or equal to 0.5
      • gDataSetsComparison

        public static double gDataSetsComparison​(long[] observed1,
                                                 long[] observed2)
                                          throws MathIllegalArgumentException

        Computes a G (Log-Likelihood Ratio) two sample test statistic for independence comparing frequency counts in observed1 and observed2. The sums of frequency counts in the two samples are not required to be the same. The formula used to compute the test statistic is

        2 * totalSum * [H(rowSums) + H(colSums) - H(k)]

        where H is the Shannon Entropy of the random variable formed by viewing the elements of the argument array as incidence counts;
        k is a matrix with rows [observed1, observed2];
        rowSums, colSums are the row/col sums of k;
        and totalSum is the overall sum of all entries in k.

        This statistic can be used to perform a G test evaluating the null hypothesis that both observed counts are independent

        Preconditions:

        • Observed counts must be non-negative.
        • Observed counts for a specific bin must not both be zero.
        • Observed counts for a specific sample must not all be 0.
        • The arrays observed1 and observed2 must have the same length and their common length must be at least 2.

        If any of the preconditions are not met, a MathIllegalArgumentException is thrown.

        Parameters:
        observed1 - array of observed frequency counts of the first data set
        observed2 - array of observed frequency counts of the second data set
        Returns:
        G-Test statistic
        Throws:
        MathIllegalArgumentException - the the lengths of the arrays do not match or their common length is less than 2
        MathIllegalArgumentException - if any entry in observed1 or observed2 is negative
        MathIllegalArgumentException - if either all counts of observed1 or observed2 are zero, or if the count at the same index is zero for both arrays.
      • rootLogLikelihoodRatio

        public static double rootLogLikelihoodRatio​(long k11,
                                                    long k12,
                                                    long k21,
                                                    long k22)
                                             throws MathIllegalArgumentException
        Calculates the root log-likelihood ratio for 2 state Datasets. See gDataSetsComparison(long[], long[] ).

        Given two events A and B, let k11 be the number of times both events occur, k12 the incidence of B without A, k21 the count of A without B, and k22 the number of times neither A nor B occurs. What is returned by this method is

        (sgn) sqrt(gValueDataSetsComparison({k11, k12}, {k21, k22})

        where sgn is -1 if k11 / (k11 + k12) < k21 / (k21 + k22));
        1 otherwise.

        Signed root LLR has two advantages over the basic LLR: a) it is positive where k11 is bigger than expected, negative where it is lower b) if there is no difference it is asymptotically normally distributed. This allows one to talk about "number of standard deviations" which is a more common frame of reference than the chi^2 distribution.

        Parameters:
        k11 - number of times the two events occurred together (AB)
        k12 - number of times the second event occurred WITHOUT the first event (notA,B)
        k21 - number of times the first event occurred WITHOUT the second event (A, notB)
        k22 - number of times something else occurred (i.e. was neither of these events (notA, notB)
        Returns:
        root log-likelihood ratio
        Throws:
        MathIllegalArgumentException
      • gTestDataSetsComparison

        public static double gTestDataSetsComparison​(long[] observed1,
                                                     long[] observed2)
                                              throws MathIllegalArgumentException,
                                                     MathIllegalStateException

        Returns the observed significance level, or p-value, associated with a G-Value (Log-Likelihood Ratio) for two sample test comparing bin frequency counts in observed1 and observed2.

        The number returned is the smallest significance level at which one can reject the null hypothesis that the observed counts conform to the same distribution.

        See gTest(double[], long[]) for details on how the p-value is computed. The degrees of of freedom used to perform the test is one less than the common length of the input observed count arrays.

        Preconditions:

        • Observed counts must be non-negative.
        • Observed counts for a specific bin must not both be zero.
        • Observed counts for a specific sample must not all be 0.
        • The arrays observed1 and observed2 must have the same length and their common length must be at least 2.

        If any of the preconditions are not met, a MathIllegalArgumentException is thrown.

        Parameters:
        observed1 - array of observed frequency counts of the first data set
        observed2 - array of observed frequency counts of the second data set
        Returns:
        p-value
        Throws:
        MathIllegalArgumentException - the the length of the arrays does not match or their common length is less than 2
        MathIllegalArgumentException - if any of the entries in observed1 or observed2 are negative
        MathIllegalArgumentException - if either all counts of observed1 or observed2 are zero, or if the count at some index is zero for both arrays
        MathIllegalStateException - if an error occurs computing the p-value.
      • gTestDataSetsComparison

        public static boolean gTestDataSetsComparison​(long[] observed1,
                                                      long[] observed2,
                                                      double alpha)
                                               throws MathIllegalArgumentException,
                                                      MathIllegalStateException

        Performs a G-Test (Log-Likelihood Ratio Test) comparing two binned data sets. The test evaluates the null hypothesis that the two lists of observed counts conform to the same frequency distribution, with significance level alpha. Returns true iff the null hypothesis can be rejected with 100 * (1 - alpha) percent confidence.

        See gDataSetsComparison(long[], long[]) for details on the formula used to compute the G (LLR) statistic used in the test and gTest(double[], long[]) for information on how the observed significance level is computed. The degrees of of freedom used to perform the test is one less than the common length of the input observed count arrays.

        Preconditions:

        • Observed counts must be non-negative.
        • Observed counts for a specific bin must not both be zero.
        • Observed counts for a specific sample must not all be 0.
        • The arrays observed1 and observed2 must have the same length and their common length must be at least 2.
        • 0 < alpha < 0.5

        If any of the preconditions are not met, a MathIllegalArgumentException is thrown.

        Parameters:
        observed1 - array of observed frequency counts of the first data set
        observed2 - array of observed frequency counts of the second data set
        alpha - significance level of the test
        Returns:
        true iff null hypothesis can be rejected with confidence 1 - alpha
        Throws:
        MathIllegalArgumentException - the the length of the arrays does not match
        MathIllegalArgumentException - if any of the entries in observed1 or observed2 are negative
        MathIllegalArgumentException - if either all counts of observed1 or observed2 are zero, or if the count at some index is zero for both arrays
        MathIllegalArgumentException - if alpha is not in the range (0, 0.5]
        MathIllegalStateException - if an error occurs performing the test
      • kolmogorovSmirnovStatistic

        public static double kolmogorovSmirnovStatistic​(RealDistribution dist,
                                                        double[] data)
                                                 throws MathIllegalArgumentException,
                                                        NullArgumentException
        Computes the one-sample Kolmogorov-Smirnov test statistic, \(D_n=\sup_x |F_n(x)-F(x)|\) where \(F\) is the distribution (cdf) function associated with distribution, \(n\) is the length of data and \(F_n\) is the empirical distribution that puts mass \(1/n\) at each of the values in data.
        Parameters:
        dist - reference distribution
        data - sample being evaluated
        Returns:
        Kolmogorov-Smirnov statistic \(D_n\)
        Throws:
        MathIllegalArgumentException - if data does not have length at least 2
        NullArgumentException - if data is null
      • kolmogorovSmirnovStatistic

        public static double kolmogorovSmirnovStatistic​(double[] x,
                                                        double[] y)
                                                 throws MathIllegalArgumentException,
                                                        NullArgumentException
        Computes the two-sample Kolmogorov-Smirnov test statistic, \(D_{n,m}=\sup_x |F_n(x)-F_m(x)|\) where \(n\) is the length of x, \(m\) is the length of y, \(F_n\) is the empirical distribution that puts mass \(1/n\) at each of the values in x and \(F_m\) is the empirical distribution of the y values.
        Parameters:
        x - first sample
        y - second sample
        Returns:
        test statistic \(D_{n,m}\) used to evaluate the null hypothesis that x and y represent samples from the same underlying distribution
        Throws:
        MathIllegalArgumentException - if either x or y does not have length at least 2
        NullArgumentException - if either x or y is null
      • exactP

        public static double exactP​(double d,
                                    int m,
                                    int n,
                                    boolean strict)
        Computes \(P(D_{n,m} > d)\) if strict is true; otherwise \(P(D_{n,m} \ge d)\), where \(D_{n,m}\) is the 2-sample Kolmogorov-Smirnov statistic. See KolmogorovSmirnovTest.kolmogorovSmirnovStatistic(double[], double[]) for the definition of \(D_{n,m}\).

        The returned probability is exact, implemented by unwinding the recursive function definitions presented in [4] from the class javadoc.

        Parameters:
        d - D-statistic value
        n - first sample size
        m - second sample size
        strict - whether or not the probability to compute is expressed as a strict inequality
        Returns:
        probability that a randomly selected m-n partition of m + n generates \(D_{n,m}\) greater than (resp. greater than or equal to) d