Quantification Metrics
This module implements some performance metrics for distribution parameterization
- class qp.metrics.metrics.Grid(grid_values, cardinality, resolution, hist_bin_edges, limits)
- cardinality
Alias for field number 1
- grid_values
Alias for field number 0
- hist_bin_edges
Alias for field number 3
- limits
Alias for field number 4
- resolution
Alias for field number 2
- qp.metrics.metrics.calculate_moment(p, N, limits, dx=0.01)[source]
Calculates a moment of a qp.Ensemble object
- Parameters:
- p: qp.Ensemble object
the collection of PDFs whose moment will be calculated
- N: int
order of the moment to be calculated
- limits: tuple of floats
endpoints of integration interval over which to calculate moments
- dx: float
resolution of integration grid
- Returns:
- M:
float value of the moment
- M:
- qp.metrics.metrics.calculate_kld(p, q, limits, dx=0.01)[source]
Calculates the Kullback-Leibler Divergence between two qp.Ensemble objects.
- Parameters:
- p: Ensemble object
probability distribution closer to the truth
- q: Ensemble object
probability distribution that approximates p
- limits: tuple of floats
endpoints of integration interval in which to calculate KLD
- dx: float
resolution of integration grid
- Returns:
- Dpq:
float the value of the Kullback-Leibler Divergence from
qtop
- Dpq:
Notes
TO DO: have this take number of points not dx!
- qp.metrics.metrics.calculate_rmse(p, q, limits, dx=0.01)[source]
Calculates the Root Mean Square Error between two qp.Ensemble objects.
- Parameters:
- p: qp.Ensemble object
probability distribution function whose distance between its truth and the approximation of
qwill be calculated.- q: qp.Ensemble object
probability distribution function whose distance between its approximation and the truth of
pwill be calculated.- limits: tuple of floats
endpoints of integration interval in which to calculate RMS
- dx: float
resolution of integration grid
- Returns:
- rms:
float the value of the RMS error between
qandp
- rms:
Notes
TO DO: change dx to N
- qp.metrics.metrics.calculate_rbpe(p, limits=(inf, inf))[source]
Calculates the risk based point estimates of a qp.Ensemble object. Algorithm as defined in 4.2 of ‘Photometric redshifts for Hyper Suprime-Cam Subaru Strategic Program Data Release 1’ (Tanaka et al. 2018).
- Parameters:
- p: qp.Ensemble object
Ensemble of PDFs to be evalutated
- limits: tuple
The limits at which to evaluate possible z_best estimates. If custom limits are not provided then all potential z value will be considered using the scipy.optimize.minimize_scalar function.
- Returns:
- rbpes:
arrayoffloats The risk based point estimates of the provided ensemble.
- rbpes:
- qp.metrics.metrics.calculate_brier(p, truth, limits, dx=0.01)[source]
This function will do the following:
Generate a Mx1 sized grid based on
limitsanddx.Produce an NxM array by evaluating the pdf for each of the N distribution objects in the Ensemble p on the grid.
Produce an NxM truth_array using the input truth and the generated grid. All values will be 0 or 1.
Create a Brier metric evaluation object
Return the result of the Brier metric calculation.
- Parameters:
- p: qp.Ensemble object
of N distributions probability distribution functions that will be gridded and compared against truth.
- truth: Nx1 sequence
the list of true values, 1 per distribution in p.
- limits: 2-tuple of floats
endpoints grid to evaluate the PDFs for the distributions in p
- dx: float
resolution of the grid Defaults to 0.01.
- Returns:
- Brier_metric:
float
- Brier_metric:
- qp.metrics.metrics.calculate_brier_for_accumulation(p, truth, limits, dx=0.01)[source]
- qp.metrics.metrics.calculate_anderson_darling(p, scipy_distribution='norm', num_samples=100, _random_state=None)[source]
This function is deprecated and will be completely removed in a later version. Please use
calculate_goodness_of_fitinstead.- Returns:
logger.warning
- qp.metrics.metrics.calculate_cramer_von_mises(p, q, num_samples=100, _random_state=None, **kwargs)[source]
This function is deprecated and will be completely removed in a later version. Please use
calculate_goodness_of_fitinstead.- Returns:
logger.warning
- qp.metrics.metrics.calculate_kolmogorov_smirnov(p, q, num_samples=100, _random_state=None)[source]
This function is deprecated and will be completely removed in a later version. Please use
calculate_goodness_of_fitinstead.- Returns:
logger.warning
- qp.metrics.metrics.calculate_outlier_rate(p, lower_limit=0.0001, upper_limit=0.9999)[source]
Fraction of outliers in each distribution
- Parameters:
- p
qp.Ensemble A collection of N distributions. This implementation expects that Ensembles are not nested.
- lower_limit
float, optional Lower bound CDF for outliers, by default 0.0001
- upper_limit
float, optional Upper bound CDF for outliers, by default 0.9999
- p
- Returns:
- [
float] 1xN array where each element is the percent of outliers for a distribution in the Ensemble.
- [
- qp.metrics.metrics.calculate_goodness_of_fit(estimate, reference, fit_metric='ks', num_samples=100, _random_state=None)[source]
This method calculates goodness of fit between the distributions in the
estimateandreferenceEnsembles using the specified fit_metric.- Parameters:
- estimate
EnsemblecontainingNdistributions Random variate samples will be drawn from this Ensemble
- reference
EnsemblecontainingNor 1 distributions The CDF of the distributions in this Ensemble are used in the goodness of fit calculation.
- fit_metric
str, optional The goodness of fit metric to use. One of [‘ad’, ‘cvm’, ‘ks’]. For clarity, ‘ad’ = Anderson-Darling, ‘cvm’ = Cramer-von Mises, and ‘ks’ = Kolmogorov-Smirnov, by default ‘ks’
- num_samples
int, optional Number of random variates to draw from each distribution in
estimate, by default 100- _random_state
_type_, optional Used for testing to create reproducible sets of random variates, by default None
- estimate
- Returns:
- output: [
float] A array of floats where each element is the result of the statistic calculation.
- output: [
- Raises:
KeyErrorIf the requested
fit_metricis not contained ingoodness_of_fit_metricsdictionary, raise a KeyError.
Notes
The calculation of the goodness of fit metrics is not symmetric. i.e.
calculate_goodness_of_fit(p, q, ...) != calculate_goodness_of_fit(q, p, ...)In the future, we should be able to do this directly from the PDFs without needing to take random variates from the
estimateEnsemble.The vectorized implementations of fit metrics are copied over (unmodified) from the developer branch of Scipy 1.10.0dev. When Scipy 1.10 is released, we can replace the copied implementation with the ones in Scipy.
This module implements metric calculations that are independent of qp.Ensembles
- qp.metrics.array_metrics.quick_anderson_ksamp(p_random_variables, q_random_variables, **kwargs)[source]
Calculate the k-sample Anderson-Darling statistic using scipy.stats.anderson_ksamp for two CDFs. For more details see: https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.anderson_ksamp.html
- Parameters:
- p_random_variables
np.array An array of random variables from the given distribution
- q_random_variables
np.array An array of random variables from the given distribution
- p_random_variables
- Returns:
- [
Resultobjects] A array of objects with attributes
statistic,critical_values, andsignificance_level.
- [
- qp.metrics.array_metrics.quick_kld(p_eval, q_eval, dx=0.01)[source]
Calculates the Kullback-Leibler Divergence between two evaluations of PDFs.
- Parameters:
- p_eval: numpy.ndarray, float
evaluations of probability distribution closer to the truth
- q_eval: numpy.ndarray, float
evaluations of probability distribution that approximates p
- dx: float
resolution of integration grid
- Returns:
- Dpq:
float the value of the Kullback-Leibler Divergence from
qtop
- Dpq:
- qp.metrics.array_metrics.quick_moment(p_eval, grid_to_N, dx)[source]
Calculates a moment of an evaluated PDF
- Parameters:
- p_eval: numpy.ndarray, float
the values of a probability distribution
- grid: numpy.ndarray, float
the grid upon which p_eval was evaluated
- dx: float
the difference between regular grid points
- N: int
order of the moment to be calculated
- Returns:
- M:
float value of the moment
- M:
- qp.metrics.array_metrics.quick_rmse(p_eval, q_eval, N)[source]
Calculates the Root Mean Square Error between two evaluations of PDFs.
- Parameters:
- p_eval: numpy.ndarray, float
evaluation of probability distribution function whose distance between its truth and the approximation of
qwill be calculated.- q_eval: numpy.ndarray, float
evaluation of probability distribution function whose distance between its approximation and the truth of
pwill be calculated.- N: int
number of points at which PDFs were evaluated
- Returns:
- rms:
float the value of the RMS error between
qandp
- rms:
- qp.metrics.array_metrics.quick_rbpe(pdf_function, integration_bounds, limits=(inf, inf))[source]
Calculates the risk based point estimate of a qp.Ensemble object with npdf == 1.
- Parameters:
- pdf_function, python function
The function should calculate the value of a pdf at a given x value
- integration_bounds, 2-tuple of floats
The integration bounds - typically (ppf(0.01), ppf(0.99)) for the given distribution
- limits, tuple of floats
The limits at which to evaluate possible z_best estimates. If custom limits are not provided then all potential z value will be considered using the scipy.optimize.minimize_scalar function.
- Returns:
- rbpe:
float The risk based point estimate of the provided ensemble.
- rbpe:
- class qp.metrics.brier.Brier(prediction, truth)[source]
Brier score based on https://en.wikipedia.org/wiki/Brier_score#Original_definition_by_Brier
- Parameters:
- prediction: NxM array, float
Predicted probability for N distributions to have a true value in one of M bins. The sum of values along each row N should be 1.
- truth: NxM array, int
True values for N distributions, where Mth bin for the true value will have value 1, all other bins will have a value of 0.
- class qp.metrics.pit.PIT(qp_ens, true_vals, eval_grid=DEFAULT_QUANTS)[source]
Probability Integral Transform
- Parameters:
- Returns:
PITobjectAn object with an Ensemble containing the PIT distribution, and a full set of PIT samples.
- property pit_samps
Returns the PIT samples. i.e.
CDF(true_vals)for each distribution in the Ensemble used to initialize the PIT object.- Returns:
np.arrayAn array of floats
- property pit
Return the PIT Ensemble object
- Returns:
qp.EnsembleAn Ensemble containing 1 qp.quant distribution.
- calculate_pit_meta_metrics()[source]
Convenience method that will calculate all of the PIT meta metrics and return them as a dictionary.
- Returns:
dictionaryThe collection of PIT statistics
- evaluate_PIT_anderson_ksamp(pit_min=0.0, pit_max=1.0)[source]
Use scipy.stats.anderson_ksamp to compute the Anderson-Darling statistic for the cdf(truth) values by comparing with a uniform distribution between 0 and 1. Up to the current version (1.9.3), scipy.stats.anderson does not support uniform distributions as reference for 1-sample test, therefore we create a uniform “distribution” and pass it as the second value in the list of parameters to the scipy implementation of k-sample Anderson-Darling. For details see: https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.anderson_ksamp.html
- Parameters:
- Returns:
arrayA array of objects with attributes
statistic,critical_values, andsignificance_level. For details see: https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.anderson_ksamp.html
- evaluate_PIT_CvM()[source]
Calculate the Cramer von Mises statistic using scipy.stats.cramervonmises using self._pit_samps compared to a uniform distribution. For more details see: https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.cramervonmises.html
- Returns:
arrayA array of objects with attributes
statisticandpvalueFor details see: https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.cramervonmises.html
- evaluate_PIT_KS()[source]
Calculate the Kolmogorov-Smirnov statistic using scipy.stats.kstest. For more details see: https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.kstest.html
- Returns:
arrayA array of objects with attributes
statisticandpvalue. For details see: https://docs.scipy.org/doc/scipy/reference/generated/scipy.stats.kstest.html
- evaluate_PIT_outlier_rate(pit_min=0.0001, pit_max=0.9999)[source]
Compute fraction of PIT outliers by evaluating the CDF of the distribution in the PIT Ensemble at
pit_minandpit_max.