_fit.py 56 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284
  1. import warnings
  2. from collections import namedtuple
  3. import numpy as np
  4. from scipy import optimize, stats
  5. from scipy._lib._util import check_random_state
  6. def _combine_bounds(name, user_bounds, shape_domain, integral):
  7. """Intersection of user-defined bounds and distribution PDF/PMF domain"""
  8. user_bounds = np.atleast_1d(user_bounds)
  9. if user_bounds[0] > user_bounds[1]:
  10. message = (f"There are no values for `{name}` on the interval "
  11. f"{list(user_bounds)}.")
  12. raise ValueError(message)
  13. bounds = (max(user_bounds[0], shape_domain[0]),
  14. min(user_bounds[1], shape_domain[1]))
  15. if integral and (np.ceil(bounds[0]) > np.floor(bounds[1])):
  16. message = (f"There are no integer values for `{name}` on the interval "
  17. f"defined by the user-provided bounds and the domain "
  18. "of the distribution.")
  19. raise ValueError(message)
  20. elif not integral and (bounds[0] > bounds[1]):
  21. message = (f"There are no values for `{name}` on the interval "
  22. f"defined by the user-provided bounds and the domain "
  23. "of the distribution.")
  24. raise ValueError(message)
  25. if not np.all(np.isfinite(bounds)):
  26. message = (f"The intersection of user-provided bounds for `{name}` "
  27. f"and the domain of the distribution is not finite. Please "
  28. f"provide finite bounds for shape `{name}` in `bounds`.")
  29. raise ValueError(message)
  30. return bounds
  31. class FitResult:
  32. r"""Result of fitting a discrete or continuous distribution to data
  33. Attributes
  34. ----------
  35. params : namedtuple
  36. A namedtuple containing the maximum likelihood estimates of the
  37. shape parameters, location, and (if applicable) scale of the
  38. distribution.
  39. success : bool or None
  40. Whether the optimizer considered the optimization to terminate
  41. successfully or not.
  42. message : str or None
  43. Any status message provided by the optimizer.
  44. """
  45. def __init__(self, dist, data, discrete, res):
  46. self._dist = dist
  47. self._data = data
  48. self.discrete = discrete
  49. self.pxf = getattr(dist, "pmf", None) or getattr(dist, "pdf", None)
  50. shape_names = [] if dist.shapes is None else dist.shapes.split(", ")
  51. if not discrete:
  52. FitParams = namedtuple('FitParams', shape_names + ['loc', 'scale'])
  53. else:
  54. FitParams = namedtuple('FitParams', shape_names + ['loc'])
  55. self.params = FitParams(*res.x)
  56. # Optimizer can report success even when nllf is infinite
  57. if res.success and not np.isfinite(self.nllf()):
  58. res.success = False
  59. res.message = ("Optimization converged to parameter values that "
  60. "are inconsistent with the data.")
  61. self.success = getattr(res, "success", None)
  62. self.message = getattr(res, "message", None)
  63. def __repr__(self):
  64. keys = ["params", "success", "message"]
  65. m = max(map(len, keys)) + 1
  66. return '\n'.join([key.rjust(m) + ': ' + repr(getattr(self, key))
  67. for key in keys if getattr(self, key) is not None])
  68. def nllf(self, params=None, data=None):
  69. """Negative log-likelihood function
  70. Evaluates the negative of the log-likelihood function of the provided
  71. data at the provided parameters.
  72. Parameters
  73. ----------
  74. params : tuple, optional
  75. The shape parameters, location, and (if applicable) scale of the
  76. distribution as a single tuple. Default is the maximum likelihood
  77. estimates (``self.params``).
  78. data : array_like, optional
  79. The data for which the log-likelihood function is to be evaluated.
  80. Default is the data to which the distribution was fit.
  81. Returns
  82. -------
  83. nllf : float
  84. The negative of the log-likelihood function.
  85. """
  86. params = params if params is not None else self.params
  87. data = data if data is not None else self._data
  88. return self._dist.nnlf(theta=params, x=data)
  89. def plot(self, ax=None, *, plot_type="hist"):
  90. """Visually compare the data against the fitted distribution.
  91. Available only if ``matplotlib`` is installed.
  92. Parameters
  93. ----------
  94. ax : matplotlib.axes.Axes
  95. Axes object to draw the plot onto, otherwise uses the current Axes.
  96. plot_type : {"hist", "qq", "pp", "cdf"}
  97. Type of plot to draw. Options include:
  98. - "hist": Superposes the PDF/PMF of the fitted distribution
  99. over a normalized histogram of the data.
  100. - "qq": Scatter plot of theoretical quantiles against the
  101. empirical quantiles. Specifically, the x-coordinates are the
  102. values of the fitted distribution PPF evaluated at the
  103. percentiles ``(np.arange(1, n) - 0.5)/n``, where ``n`` is the
  104. number of data points, and the y-coordinates are the sorted
  105. data points.
  106. - "pp": Scatter plot of theoretical percentiles against the
  107. observed percentiles. Specifically, the x-coordinates are the
  108. percentiles ``(np.arange(1, n) - 0.5)/n``, where ``n`` is
  109. the number of data points, and the y-coordinates are the values
  110. of the fitted distribution CDF evaluated at the sorted
  111. data points.
  112. - "cdf": Superposes the CDF of the fitted distribution over the
  113. empirical CDF. Specifically, the x-coordinates of the empirical
  114. CDF are the sorted data points, and the y-coordinates are the
  115. percentiles ``(np.arange(1, n) - 0.5)/n``, where ``n`` is
  116. the number of data points.
  117. Returns
  118. -------
  119. ax : matplotlib.axes.Axes
  120. The matplotlib Axes object on which the plot was drawn.
  121. """
  122. try:
  123. import matplotlib # noqa
  124. except ModuleNotFoundError as exc:
  125. message = "matplotlib must be installed to use method `plot`."
  126. raise ModuleNotFoundError(message) from exc
  127. plots = {'histogram': self._hist_plot, 'qq': self._qq_plot,
  128. 'pp': self._pp_plot, 'cdf': self._cdf_plot,
  129. 'hist': self._hist_plot}
  130. if plot_type.lower() not in plots:
  131. message = f"`plot_type` must be one of {set(plots.keys())}"
  132. raise ValueError(message)
  133. plot = plots[plot_type.lower()]
  134. if ax is None:
  135. import matplotlib.pyplot as plt
  136. ax = plt.gca()
  137. fit_params = np.atleast_1d(self.params)
  138. return plot(ax=ax, fit_params=fit_params)
  139. def _hist_plot(self, ax, fit_params):
  140. from matplotlib.ticker import MaxNLocator
  141. support = self._dist.support(*fit_params)
  142. lb = support[0] if np.isfinite(support[0]) else min(self._data)
  143. ub = support[1] if np.isfinite(support[1]) else max(self._data)
  144. pxf = "PMF" if self.discrete else "PDF"
  145. if self.discrete:
  146. x = np.arange(lb, ub + 2)
  147. y = self.pxf(x, *fit_params)
  148. ax.vlines(x[:-1], 0, y[:-1], label='Fitted Distribution PMF',
  149. color='C0')
  150. options = dict(density=True, bins=x, align='left', color='C1')
  151. ax.xaxis.set_major_locator(MaxNLocator(integer=True))
  152. ax.set_xlabel('k')
  153. ax.set_ylabel('PMF')
  154. else:
  155. x = np.linspace(lb, ub, 200)
  156. y = self.pxf(x, *fit_params)
  157. ax.plot(x, y, '--', label='Fitted Distribution PDF', color='C0')
  158. options = dict(density=True, bins=50, align='mid', color='C1')
  159. ax.set_xlabel('x')
  160. ax.set_ylabel('PDF')
  161. if len(self._data) > 50 or self.discrete:
  162. ax.hist(self._data, label="Histogram of Data", **options)
  163. else:
  164. ax.plot(self._data, np.zeros_like(self._data), "*",
  165. label='Data', color='C1')
  166. ax.set_title(rf"Fitted $\tt {self._dist.name}$ {pxf} and Histogram")
  167. ax.legend(*ax.get_legend_handles_labels())
  168. return ax
  169. def _qp_plot(self, ax, fit_params, qq):
  170. data = np.sort(self._data)
  171. ps = self._plotting_positions(len(self._data))
  172. if qq:
  173. qp = "Quantiles"
  174. plot_type = 'Q-Q'
  175. x = self._dist.ppf(ps, *fit_params)
  176. y = data
  177. else:
  178. qp = "Percentiles"
  179. plot_type = 'P-P'
  180. x = ps
  181. y = self._dist.cdf(data, *fit_params)
  182. ax.plot(x, y, '.', label=f'Fitted Distribution {plot_type}',
  183. color='C0', zorder=1)
  184. xlim = ax.get_xlim()
  185. ylim = ax.get_ylim()
  186. lim = [min(xlim[0], ylim[0]), max(xlim[1], ylim[1])]
  187. if not qq:
  188. lim = max(lim[0], 0), min(lim[1], 1)
  189. if self.discrete and qq:
  190. q_min, q_max = int(lim[0]), int(lim[1]+1)
  191. q_ideal = np.arange(q_min, q_max)
  192. # q_ideal = np.unique(self._dist.ppf(ps, *fit_params))
  193. ax.plot(q_ideal, q_ideal, 'o', label='Reference', color='k',
  194. alpha=0.25, markerfacecolor='none', clip_on=True)
  195. elif self.discrete and not qq:
  196. # The intent of this is to match the plot that would be produced
  197. # if x were continuous on [0, 1] and y were cdf(ppf(x)).
  198. # It can be approximated by letting x = np.linspace(0, 1, 1000),
  199. # but this might not look great when zooming in. The vertical
  200. # portions are included to indicate where the transition occurs
  201. # where the data completely obscures the horizontal portions.
  202. p_min, p_max = lim
  203. a, b = self._dist.support(*fit_params)
  204. p_min = max(p_min, 0 if np.isfinite(a) else 1e-3)
  205. p_max = min(p_max, 1 if np.isfinite(b) else 1-1e-3)
  206. q_min, q_max = self._dist.ppf([p_min, p_max], *fit_params)
  207. qs = np.arange(q_min-1, q_max+1)
  208. ps = self._dist.cdf(qs, *fit_params)
  209. ax.step(ps, ps, '-', label='Reference', color='k', alpha=0.25,
  210. clip_on=True)
  211. else:
  212. ax.plot(lim, lim, '-', label='Reference', color='k', alpha=0.25,
  213. clip_on=True)
  214. ax.set_xlim(lim)
  215. ax.set_ylim(lim)
  216. ax.set_xlabel(rf"Fitted $\tt {self._dist.name}$ Theoretical {qp}")
  217. ax.set_ylabel(f"Data {qp}")
  218. ax.set_title(rf"Fitted $\tt {self._dist.name}$ {plot_type} Plot")
  219. ax.legend(*ax.get_legend_handles_labels())
  220. ax.set_aspect('equal')
  221. return ax
  222. def _qq_plot(self, **kwargs):
  223. return self._qp_plot(qq=True, **kwargs)
  224. def _pp_plot(self, **kwargs):
  225. return self._qp_plot(qq=False, **kwargs)
  226. def _plotting_positions(self, n, a=.5):
  227. # See https://en.wikipedia.org/wiki/Q%E2%80%93Q_plot#Plotting_positions
  228. k = np.arange(1, n+1)
  229. return (k-a) / (n + 1 - 2*a)
  230. def _cdf_plot(self, ax, fit_params):
  231. data = np.sort(self._data)
  232. ecdf = self._plotting_positions(len(self._data))
  233. ls = '--' if len(np.unique(data)) < 30 else '.'
  234. xlabel = 'k' if self.discrete else 'x'
  235. ax.step(data, ecdf, ls, label='Empirical CDF', color='C1', zorder=0)
  236. xlim = ax.get_xlim()
  237. q = np.linspace(*xlim, 300)
  238. tcdf = self._dist.cdf(q, *fit_params)
  239. ax.plot(q, tcdf, label='Fitted Distribution CDF', color='C0', zorder=1)
  240. ax.set_xlim(xlim)
  241. ax.set_ylim(0, 1)
  242. ax.set_xlabel(xlabel)
  243. ax.set_ylabel("CDF")
  244. ax.set_title(rf"Fitted $\tt {self._dist.name}$ and Empirical CDF")
  245. handles, labels = ax.get_legend_handles_labels()
  246. ax.legend(handles[::-1], labels[::-1])
  247. return ax
  248. def fit(dist, data, bounds=None, *, guess=None, method='mle',
  249. optimizer=optimize.differential_evolution):
  250. r"""Fit a discrete or continuous distribution to data
  251. Given a distribution, data, and bounds on the parameters of the
  252. distribution, return maximum likelihood estimates of the parameters.
  253. Parameters
  254. ----------
  255. dist : `scipy.stats.rv_continuous` or `scipy.stats.rv_discrete`
  256. The object representing the distribution to be fit to the data.
  257. data : 1D array_like
  258. The data to which the distribution is to be fit. If the data contain
  259. any of ``np.nan``, ``np.inf``, or -``np.inf``, the fit method will
  260. raise a ``ValueError``.
  261. bounds : dict or sequence of tuples, optional
  262. If a dictionary, each key is the name of a parameter of the
  263. distribution, and the corresponding value is a tuple containing the
  264. lower and upper bound on that parameter. If the distribution is
  265. defined only for a finite range of values of that parameter, no entry
  266. for that parameter is required; e.g., some distributions have
  267. parameters which must be on the interval [0, 1]. Bounds for parameters
  268. location (``loc``) and scale (``scale``) are optional; by default,
  269. they are fixed to 0 and 1, respectively.
  270. If a sequence, element *i* is a tuple containing the lower and upper
  271. bound on the *i*\ th parameter of the distribution. In this case,
  272. bounds for *all* distribution shape parameters must be provided.
  273. Optionally, bounds for location and scale may follow the
  274. distribution shape parameters.
  275. If a shape is to be held fixed (e.g. if it is known), the
  276. lower and upper bounds may be equal. If a user-provided lower or upper
  277. bound is beyond a bound of the domain for which the distribution is
  278. defined, the bound of the distribution's domain will replace the
  279. user-provided value. Similarly, parameters which must be integral
  280. will be constrained to integral values within the user-provided bounds.
  281. guess : dict or array_like, optional
  282. If a dictionary, each key is the name of a parameter of the
  283. distribution, and the corresponding value is a guess for the value
  284. of the parameter.
  285. If a sequence, element *i* is a guess for the *i*\ th parameter of the
  286. distribution. In this case, guesses for *all* distribution shape
  287. parameters must be provided.
  288. If `guess` is not provided, guesses for the decision variables will
  289. not be passed to the optimizer. If `guess` is provided, guesses for
  290. any missing parameters will be set at the mean of the lower and
  291. upper bounds. Guesses for parameters which must be integral will be
  292. rounded to integral values, and guesses that lie outside the
  293. intersection of the user-provided bounds and the domain of the
  294. distribution will be clipped.
  295. method : {'mle', 'mse'}
  296. With ``method="mle"`` (default), the fit is computed by minimizing
  297. the negative log-likelihood function. A large, finite penalty
  298. (rather than infinite negative log-likelihood) is applied for
  299. observations beyond the support of the distribution.
  300. With ``method="mse"``, the fit is computed by minimizing
  301. the negative log-product spacing function. The same penalty is applied
  302. for observations beyond the support. We follow the approach of [1]_,
  303. which is generalized for samples with repeated observations.
  304. optimizer : callable, optional
  305. `optimizer` is a callable that accepts the following positional
  306. argument.
  307. fun : callable
  308. The objective function to be optimized. `fun` accepts one argument
  309. ``x``, candidate shape parameters of the distribution, and returns
  310. the objective function value given ``x``, `dist`, and the provided
  311. `data`.
  312. The job of `optimizer` is to find values of the decision variables
  313. that minimizes `fun`.
  314. `optimizer` must also accept the following keyword argument.
  315. bounds : sequence of tuples
  316. The bounds on values of the decision variables; each element will
  317. be a tuple containing the lower and upper bound on a decision
  318. variable.
  319. If `guess` is provided, `optimizer` must also accept the following
  320. keyword argument.
  321. x0 : array_like
  322. The guesses for each decision variable.
  323. If the distribution has any shape parameters that must be integral or
  324. if the distribution is discrete and the location parameter is not
  325. fixed, `optimizer` must also accept the following keyword argument.
  326. integrality : array_like of bools
  327. For each decision variable, True if the decision variable
  328. must be constrained to integer values and False if the decision
  329. variable is continuous.
  330. `optimizer` must return an object, such as an instance of
  331. `scipy.optimize.OptimizeResult`, which holds the optimal values of
  332. the decision variables in an attribute ``x``. If attributes
  333. ``fun``, ``status``, or ``message`` are provided, they will be
  334. included in the result object returned by `fit`.
  335. Returns
  336. -------
  337. result : `~scipy.stats._result_classes.FitResult`
  338. An object with the following fields.
  339. params : namedtuple
  340. A namedtuple containing the maximum likelihood estimates of the
  341. shape parameters, location, and (if applicable) scale of the
  342. distribution.
  343. success : bool or None
  344. Whether the optimizer considered the optimization to terminate
  345. successfully or not.
  346. message : str or None
  347. Any status message provided by the optimizer.
  348. The object has the following method:
  349. nllf(params=None, data=None)
  350. By default, the negative log-likehood function at the fitted
  351. `params` for the given `data`. Accepts a tuple containing
  352. alternative shapes, location, and scale of the distribution and
  353. an array of alternative data.
  354. plot(ax=None)
  355. Superposes the PDF/PMF of the fitted distribution over a normalized
  356. histogram of the data.
  357. See Also
  358. --------
  359. rv_continuous, rv_discrete
  360. Notes
  361. -----
  362. Optimization is more likely to converge to the maximum likelihood estimate
  363. when the user provides tight bounds containing the maximum likelihood
  364. estimate. For example, when fitting a binomial distribution to data, the
  365. number of experiments underlying each sample may be known, in which case
  366. the corresponding shape parameter ``n`` can be fixed.
  367. References
  368. ----------
  369. .. [1] Shao, Yongzhao, and Marjorie G. Hahn. "Maximum product of spacings
  370. method: a unified formulation with illustration of strong
  371. consistency." Illinois Journal of Mathematics 43.3 (1999): 489-499.
  372. Examples
  373. --------
  374. Suppose we wish to fit a distribution to the following data.
  375. >>> import numpy as np
  376. >>> from scipy import stats
  377. >>> rng = np.random.default_rng()
  378. >>> dist = stats.nbinom
  379. >>> shapes = (5, 0.5)
  380. >>> data = dist.rvs(*shapes, size=1000, random_state=rng)
  381. Suppose we do not know how the data were generated, but we suspect that
  382. it follows a negative binomial distribution with parameters *n* and *p*\.
  383. (See `scipy.stats.nbinom`.) We believe that the parameter *n* was fewer
  384. than 30, and we know that the parameter *p* must lie on the interval
  385. [0, 1]. We record this information in a variable `bounds` and pass
  386. this information to `fit`.
  387. >>> bounds = [(0, 30), (0, 1)]
  388. >>> res = stats.fit(dist, data, bounds)
  389. `fit` searches within the user-specified `bounds` for the
  390. values that best match the data (in the sense of maximum likelihood
  391. estimation). In this case, it found shape values similar to those
  392. from which the data were actually generated.
  393. >>> res.params
  394. FitParams(n=5.0, p=0.5028157644634368, loc=0.0) # may vary
  395. We can visualize the results by superposing the probability mass function
  396. of the distribution (with the shapes fit to the data) over a normalized
  397. histogram of the data.
  398. >>> import matplotlib.pyplot as plt # matplotlib must be installed to plot
  399. >>> res.plot()
  400. >>> plt.show()
  401. Note that the estimate for *n* was exactly integral; this is because
  402. the domain of the `nbinom` PMF includes only integral *n*, and the `nbinom`
  403. object "knows" that. `nbinom` also knows that the shape *p* must be a
  404. value between 0 and 1. In such a case - when the domain of the distribution
  405. with respect to a parameter is finite - we are not required to specify
  406. bounds for the parameter.
  407. >>> bounds = {'n': (0, 30)} # omit parameter p using a `dict`
  408. >>> res2 = stats.fit(dist, data, bounds)
  409. >>> res2.params
  410. FitParams(n=5.0, p=0.5016492009232932, loc=0.0) # may vary
  411. If we wish to force the distribution to be fit with *n* fixed at 6, we can
  412. set both the lower and upper bounds on *n* to 6. Note, however, that the
  413. value of the objective function being optimized is typically worse (higher)
  414. in this case.
  415. >>> bounds = {'n': (6, 6)} # fix parameter `n`
  416. >>> res3 = stats.fit(dist, data, bounds)
  417. >>> res3.params
  418. FitParams(n=6.0, p=0.5486556076755706, loc=0.0) # may vary
  419. >>> res3.nllf() > res.nllf()
  420. True # may vary
  421. Note that the numerical results of the previous examples are typical, but
  422. they may vary because the default optimizer used by `fit`,
  423. `scipy.optimize.differential_evolution`, is stochastic. However, we can
  424. customize the settings used by the optimizer to ensure reproducibility -
  425. or even use a different optimizer entirely - using the `optimizer`
  426. parameter.
  427. >>> from scipy.optimize import differential_evolution
  428. >>> rng = np.random.default_rng(767585560716548)
  429. >>> def optimizer(fun, bounds, *, integrality):
  430. ... return differential_evolution(fun, bounds, strategy='best2bin',
  431. ... seed=rng, integrality=integrality)
  432. >>> bounds = [(0, 30), (0, 1)]
  433. >>> res4 = stats.fit(dist, data, bounds, optimizer=optimizer)
  434. >>> res4.params
  435. FitParams(n=5.0, p=0.5015183149259951, loc=0.0)
  436. """
  437. # --- Input Validation / Standardization --- #
  438. user_bounds = bounds
  439. user_guess = guess
  440. # distribution input validation and information collection
  441. if hasattr(dist, "pdf"): # can't use isinstance for types
  442. default_bounds = {'loc': (0, 0), 'scale': (1, 1)}
  443. discrete = False
  444. elif hasattr(dist, "pmf"):
  445. default_bounds = {'loc': (0, 0)}
  446. discrete = True
  447. else:
  448. message = ("`dist` must be an instance of `rv_continuous` "
  449. "or `rv_discrete.`")
  450. raise ValueError(message)
  451. try:
  452. param_info = dist._param_info()
  453. except AttributeError as e:
  454. message = (f"Distribution `{dist.name}` is not yet supported by "
  455. "`scipy.stats.fit` because shape information has "
  456. "not been defined.")
  457. raise ValueError(message) from e
  458. # data input validation
  459. data = np.asarray(data)
  460. if data.ndim != 1:
  461. message = "`data` must be exactly one-dimensional."
  462. raise ValueError(message)
  463. if not (np.issubdtype(data.dtype, np.number)
  464. and np.all(np.isfinite(data))):
  465. message = "All elements of `data` must be finite numbers."
  466. raise ValueError(message)
  467. # bounds input validation and information collection
  468. n_params = len(param_info)
  469. n_shapes = n_params - (1 if discrete else 2)
  470. param_list = [param.name for param in param_info]
  471. param_names = ", ".join(param_list)
  472. shape_names = ", ".join(param_list[:n_shapes])
  473. if user_bounds is None:
  474. user_bounds = {}
  475. if isinstance(user_bounds, dict):
  476. default_bounds.update(user_bounds)
  477. user_bounds = default_bounds
  478. user_bounds_array = np.empty((n_params, 2))
  479. for i in range(n_params):
  480. param_name = param_info[i].name
  481. user_bound = user_bounds.pop(param_name, None)
  482. if user_bound is None:
  483. user_bound = param_info[i].domain
  484. user_bounds_array[i] = user_bound
  485. if user_bounds:
  486. message = ("Bounds provided for the following unrecognized "
  487. f"parameters will be ignored: {set(user_bounds)}")
  488. warnings.warn(message, RuntimeWarning, stacklevel=2)
  489. else:
  490. try:
  491. user_bounds = np.asarray(user_bounds, dtype=float)
  492. if user_bounds.size == 0:
  493. user_bounds = np.empty((0, 2))
  494. except ValueError as e:
  495. message = ("Each element of a `bounds` sequence must be a tuple "
  496. "containing two elements: the lower and upper bound of "
  497. "a distribution parameter.")
  498. raise ValueError(message) from e
  499. if (user_bounds.ndim != 2 or user_bounds.shape[1] != 2):
  500. message = ("Each element of `bounds` must be a tuple specifying "
  501. "the lower and upper bounds of a shape parameter")
  502. raise ValueError(message)
  503. if user_bounds.shape[0] < n_shapes:
  504. message = (f"A `bounds` sequence must contain at least {n_shapes} "
  505. "elements: tuples specifying the lower and upper "
  506. f"bounds of all shape parameters {shape_names}.")
  507. raise ValueError(message)
  508. if user_bounds.shape[0] > n_params:
  509. message = ("A `bounds` sequence may not contain more than "
  510. f"{n_params} elements: tuples specifying the lower and "
  511. "upper bounds of distribution parameters "
  512. f"{param_names}.")
  513. raise ValueError(message)
  514. user_bounds_array = np.empty((n_params, 2))
  515. user_bounds_array[n_shapes:] = list(default_bounds.values())
  516. user_bounds_array[:len(user_bounds)] = user_bounds
  517. user_bounds = user_bounds_array
  518. validated_bounds = []
  519. for i in range(n_params):
  520. name = param_info[i].name
  521. user_bound = user_bounds_array[i]
  522. param_domain = param_info[i].domain
  523. integral = param_info[i].integrality
  524. combined = _combine_bounds(name, user_bound, param_domain, integral)
  525. validated_bounds.append(combined)
  526. bounds = np.asarray(validated_bounds)
  527. integrality = [param.integrality for param in param_info]
  528. # guess input validation
  529. if user_guess is None:
  530. guess_array = None
  531. elif isinstance(user_guess, dict):
  532. default_guess = {param.name: np.mean(bound)
  533. for param, bound in zip(param_info, bounds)}
  534. unrecognized = set(user_guess) - set(default_guess)
  535. if unrecognized:
  536. message = ("Guesses provided for the following unrecognized "
  537. f"parameters will be ignored: {unrecognized}")
  538. warnings.warn(message, RuntimeWarning, stacklevel=2)
  539. default_guess.update(user_guess)
  540. message = ("Each element of `guess` must be a scalar "
  541. "guess for a distribution parameter.")
  542. try:
  543. guess_array = np.asarray([default_guess[param.name]
  544. for param in param_info], dtype=float)
  545. except ValueError as e:
  546. raise ValueError(message) from e
  547. else:
  548. message = ("Each element of `guess` must be a scalar "
  549. "guess for a distribution parameter.")
  550. try:
  551. user_guess = np.asarray(user_guess, dtype=float)
  552. except ValueError as e:
  553. raise ValueError(message) from e
  554. if user_guess.ndim != 1:
  555. raise ValueError(message)
  556. if user_guess.shape[0] < n_shapes:
  557. message = (f"A `guess` sequence must contain at least {n_shapes} "
  558. "elements: scalar guesses for the distribution shape "
  559. f"parameters {shape_names}.")
  560. raise ValueError(message)
  561. if user_guess.shape[0] > n_params:
  562. message = ("A `guess` sequence may not contain more than "
  563. f"{n_params} elements: scalar guesses for the "
  564. f"distribution parameters {param_names}.")
  565. raise ValueError(message)
  566. guess_array = np.mean(bounds, axis=1)
  567. guess_array[:len(user_guess)] = user_guess
  568. if guess_array is not None:
  569. guess_rounded = guess_array.copy()
  570. guess_rounded[integrality] = np.round(guess_rounded[integrality])
  571. rounded = np.where(guess_rounded != guess_array)[0]
  572. for i in rounded:
  573. message = (f"Guess for parameter `{param_info[i].name}` "
  574. f"rounded from {guess_array[i]} to {guess_rounded[i]}.")
  575. warnings.warn(message, RuntimeWarning, stacklevel=2)
  576. guess_clipped = np.clip(guess_rounded, bounds[:, 0], bounds[:, 1])
  577. clipped = np.where(guess_clipped != guess_rounded)[0]
  578. for i in clipped:
  579. message = (f"Guess for parameter `{param_info[i].name}` "
  580. f"clipped from {guess_rounded[i]} to "
  581. f"{guess_clipped[i]}.")
  582. warnings.warn(message, RuntimeWarning, stacklevel=2)
  583. guess = guess_clipped
  584. else:
  585. guess = None
  586. # --- Fitting --- #
  587. def nllf(free_params, data=data): # bind data NOW
  588. with np.errstate(invalid='ignore', divide='ignore'):
  589. return dist._penalized_nnlf(free_params, data)
  590. def nlpsf(free_params, data=data): # bind data NOW
  591. with np.errstate(invalid='ignore', divide='ignore'):
  592. return dist._penalized_nlpsf(free_params, data)
  593. methods = {'mle': nllf, 'mse': nlpsf}
  594. objective = methods[method.lower()]
  595. with np.errstate(invalid='ignore', divide='ignore'):
  596. kwds = {}
  597. if bounds is not None:
  598. kwds['bounds'] = bounds
  599. if np.any(integrality):
  600. kwds['integrality'] = integrality
  601. if guess is not None:
  602. kwds['x0'] = guess
  603. res = optimizer(objective, **kwds)
  604. return FitResult(dist, data, discrete, res)
  605. GoodnessOfFitResult = namedtuple('GoodnessOfFitResult',
  606. ('fit_result', 'statistic', 'pvalue',
  607. 'null_distribution'))
  608. def goodness_of_fit(dist, data, *, known_params=None, fit_params=None,
  609. guessed_params=None, statistic='ad', n_mc_samples=9999,
  610. random_state=None):
  611. r"""
  612. Perform a goodness of fit test comparing data to a distribution family.
  613. Given a distribution family and data, perform a test of the null hypothesis
  614. that the data were drawn from a distribution in that family. Any known
  615. parameters of the distribution may be specified. Remaining parameters of
  616. the distribution will be fit to the data, and the p-value of the test
  617. is computed accordingly. Several statistics for comparing the distribution
  618. to data are available.
  619. Parameters
  620. ----------
  621. dist : `scipy.stats.rv_continuous`
  622. The object representing the distribution family under the null
  623. hypothesis.
  624. data : 1D array_like
  625. Finite, uncensored data to be tested.
  626. known_params : dict, optional
  627. A dictionary containing name-value pairs of known distribution
  628. parameters. Monte Carlo samples are randomly drawn from the
  629. null-hypothesized distribution with these values of the parameters.
  630. Before the statistic is evaluated for each Monte Carlo sample, only
  631. remaining unknown parameters of the null-hypothesized distribution
  632. family are fit to the samples; the known parameters are held fixed.
  633. If all parameters of the distribution family are known, then the step
  634. of fitting the distribution family to each sample is omitted.
  635. fit_params : dict, optional
  636. A dictionary containing name-value pairs of distribution parameters
  637. that have already been fit to the data, e.g. using `scipy.stats.fit`
  638. or the ``fit`` method of `dist`. Monte Carlo samples are drawn from the
  639. null-hypothesized distribution with these specified values of the
  640. parameter. On those Monte Carlo samples, however, these and all other
  641. unknown parameters of the null-hypothesized distribution family are
  642. fit before the statistic is evaluated.
  643. guessed_params : dict, optional
  644. A dictionary containing name-value pairs of distribution parameters
  645. which have been guessed. These parameters are always considered as
  646. free parameters and are fit both to the provided `data` as well as
  647. to the Monte Carlo samples drawn from the null-hypothesized
  648. distribution. The purpose of these `guessed_params` is to be used as
  649. initial values for the numerical fitting procedure.
  650. statistic : {"ad", "ks", "cvm"}, optional
  651. The statistic used to compare data to a distribution after fitting
  652. unknown parameters of the distribution family to the data. The
  653. Anderson-Darling ("ad"), Kolmogorov-Smirnov ("ks"), and
  654. Cramer-von Mises ("cvm") statistics are available [1]_.
  655. n_mc_samples : int, default: 9999
  656. The number of Monte Carlo samples drawn from the null hypothesized
  657. distribution to form the null distribution of the statistic. The
  658. sample size of each is the same as the given `data`.
  659. random_state : {None, int, `numpy.random.Generator`,
  660. `numpy.random.RandomState`}, optional
  661. Pseudorandom number generator state used to generate the Monte Carlo
  662. samples.
  663. If `random_state` is ``None`` (default), the
  664. `numpy.random.RandomState` singleton is used.
  665. If `random_state` is an int, a new ``RandomState`` instance is used,
  666. seeded with `random_state`.
  667. If `random_state` is already a ``Generator`` or ``RandomState``
  668. instance, then the provided instance is used.
  669. Returns
  670. -------
  671. res : GoodnessOfFitResult
  672. An object with the following attributes.
  673. fit_result : `~scipy.stats._result_classes.FitResult`
  674. An object representing the fit of the provided `dist` to `data`.
  675. This object includes the values of distribution family parameters
  676. that fully define the null-hypothesized distribution, that is,
  677. the distribution from which Monte Carlo samples are drawn.
  678. statistic : float
  679. The value of the statistic comparing provided `data` to the
  680. null-hypothesized distribution.
  681. pvalue : float
  682. The proportion of elements in the null distribution with
  683. statistic values at least as extreme as the statistic value of the
  684. provided `data`.
  685. null_distribution : ndarray
  686. The value of the statistic for each Monte Carlo sample
  687. drawn from the null-hypothesized distribution.
  688. Notes
  689. -----
  690. This is a generalized Monte Carlo goodness-of-fit procedure, special cases
  691. of which correspond with various Anderson-Darling tests, Lilliefors' test,
  692. etc. The test is described in [2]_, [3]_, and [4]_ as a parametric
  693. bootstrap test. This is a Monte Carlo test in which parameters that
  694. specify the distribution from which samples are drawn have been estimated
  695. from the data. We describe the test using "Monte Carlo" rather than
  696. "parametric bootstrap" throughout to avoid confusion with the more familiar
  697. nonparametric bootstrap, and describe how the test is performed below.
  698. *Traditional goodness of fit tests*
  699. Traditionally, critical values corresponding with a fixed set of
  700. significance levels are pre-calculated using Monte Carlo methods. Users
  701. perform the test by calculating the value of the test statistic only for
  702. their observed `data` and comparing this value to tabulated critical
  703. values. This practice is not very flexible, as tables are not available for
  704. all distributions and combinations of known and unknown parameter values.
  705. Also, results can be inaccurate when critical values are interpolated from
  706. limited tabulated data to correspond with the user's sample size and
  707. fitted parameter values. To overcome these shortcomings, this function
  708. allows the user to perform the Monte Carlo trials adapted to their
  709. particular data.
  710. *Algorithmic overview*
  711. In brief, this routine executes the following steps:
  712. 1. Fit unknown parameters to the given `data`, thereby forming the
  713. "null-hypothesized" distribution, and compute the statistic of
  714. this pair of data and distribution.
  715. 2. Draw random samples from this null-hypothesized distribution.
  716. 3. Fit the unknown parameters to each random sample.
  717. 4. Calculate the statistic between each sample and the distribution that
  718. has been fit to the sample.
  719. 5. Compare the value of the statistic corresponding with `data` from (1)
  720. against the values of the statistic corresponding with the random
  721. samples from (4). The p-value is the proportion of samples with a
  722. statistic value greater than or equal to the statistic of the observed
  723. data.
  724. In more detail, the steps are as follows.
  725. First, any unknown parameters of the distribution family specified by
  726. `dist` are fit to the provided `data` using maximum likelihood estimation.
  727. (One exception is the normal distribution with unknown location and scale:
  728. we use the bias-corrected standard deviation ``np.std(data, ddof=1)`` for
  729. the scale as recommended in [1]_.)
  730. These values of the parameters specify a particular member of the
  731. distribution family referred to as the "null-hypothesized distribution",
  732. that is, the distribution from which the data were sampled under the null
  733. hypothesis. The `statistic`, which compares data to a distribution, is
  734. computed between `data` and the null-hypothesized distribution.
  735. Next, many (specifically `n_mc_samples`) new samples, each containing the
  736. same number of observations as `data`, are drawn from the
  737. null-hypothesized distribution. All unknown parameters of the distribution
  738. family `dist` are fit to *each resample*, and the `statistic` is computed
  739. between each sample and its corresponding fitted distribution. These
  740. values of the statistic form the Monte Carlo null distribution (not to be
  741. confused with the "null-hypothesized distribution" above).
  742. The p-value of the test is the proportion of statistic values in the Monte
  743. Carlo null distribution that are at least as extreme as the statistic value
  744. of the provided `data`. More precisely, the p-value is given by
  745. .. math::
  746. p = \frac{b + 1}
  747. {m + 1}
  748. where :math:`b` is the number of statistic values in the Monte Carlo null
  749. distribution that are greater than or equal to the the statistic value
  750. calculated for `data`, and :math:`m` is the number of elements in the
  751. Monte Carlo null distribution (`n_mc_samples`). The addition of :math:`1`
  752. to the numerator and denominator can be thought of as including the
  753. value of the statistic corresponding with `data` in the null distribution,
  754. but a more formal explanation is given in [5]_.
  755. *Limitations*
  756. The test can be very slow for some distribution families because unknown
  757. parameters of the distribution family must be fit to each of the Monte
  758. Carlo samples, and for most distributions in SciPy, distribution fitting
  759. performed via numerical optimization.
  760. *Anti-Pattern*
  761. For this reason, it may be tempting
  762. to treat parameters of the distribution pre-fit to `data` (by the user)
  763. as though they were `known_params`, as specification of all parameters of
  764. the distribution precludes the need to fit the distribution to each Monte
  765. Carlo sample. (This is essentially how the original Kilmogorov-Smirnov
  766. test is performed.) Although such a test can provide evidence against the
  767. null hypothesis, the test is conservative in the sense that small p-values
  768. will tend to (greatly) *overestimate* the probability of making a type I
  769. error (that is, rejecting the null hypothesis although it is true), and the
  770. power of the test is low (that is, it is less likely to reject the null
  771. hypothesis even when the null hypothesis is false).
  772. This is because the Monte Carlo samples are less likely to agree with the
  773. null-hypothesized distribution as well as `data`. This tends to increase
  774. the values of the statistic recorded in the null distribution, so that a
  775. larger number of them exceed the value of statistic for `data`, thereby
  776. inflating the p-value.
  777. References
  778. ----------
  779. .. [1] M. A. Stephens (1974). "EDF Statistics for Goodness of Fit and
  780. Some Comparisons." Journal of the American Statistical Association,
  781. Vol. 69, pp. 730-737.
  782. .. [2] W. Stute, W. G. Manteiga, and M. P. Quindimil (1993).
  783. "Bootstrap based goodness-of-fit-tests." Metrika 40.1: 243-256.
  784. .. [3] C. Genest, & B Rémillard. (2008). "Validity of the parametric
  785. bootstrap for goodness-of-fit testing in semiparametric models."
  786. Annales de l'IHP Probabilités et statistiques. Vol. 44. No. 6.
  787. .. [4] I. Kojadinovic and J. Yan (2012). "Goodness-of-fit testing based on
  788. a weighted bootstrap: A fast large-sample alternative to the
  789. parametric bootstrap." Canadian Journal of Statistics 40.3: 480-500.
  790. .. [5] B. Phipson and G. K. Smyth (2010). "Permutation P-values Should
  791. Never Be Zero: Calculating Exact P-values When Permutations Are
  792. Randomly Drawn." Statistical Applications in Genetics and Molecular
  793. Biology 9.1.
  794. .. [6] H. W. Lilliefors (1967). "On the Kolmogorov-Smirnov test for
  795. normality with mean and variance unknown." Journal of the American
  796. statistical Association 62.318: 399-402.
  797. Examples
  798. --------
  799. A well-known test of the null hypothesis that data were drawn from a
  800. given distribution is the Kolmogorov-Smirnov (KS) test, available in SciPy
  801. as `scipy.stats.ks_1samp`. Suppose we wish to test whether the following
  802. data:
  803. >>> import numpy as np
  804. >>> from scipy import stats
  805. >>> rng = np.random.default_rng()
  806. >>> x = stats.uniform.rvs(size=75, random_state=rng)
  807. were sampled from a normal distribution. To perform a KS test, the
  808. empirical distribution function of the observed data will be compared
  809. against the (theoretical) cumulative distribution function of a normal
  810. distribution. Of course, to do this, the normal distribution under the null
  811. hypothesis must be fully specified. This is commonly done by first fitting
  812. the ``loc`` and ``scale`` parameters of the distribution to the observed
  813. data, then performing the test.
  814. >>> loc, scale = np.mean(x), np.std(x, ddof=1)
  815. >>> cdf = stats.norm(loc, scale).cdf
  816. >>> stats.ks_1samp(x, cdf)
  817. KstestResult(statistic=0.1119257570456813, pvalue=0.2827756409939257)
  818. An advantage of the KS-test is that the p-value - the probability of
  819. obtaining a value of the test statistic under the null hypothesis as
  820. extreme as the value obtained from the observed data - can be calculated
  821. exactly and efficiently. `goodness_of_fit` can only approximate these
  822. results.
  823. >>> known_params = {'loc': loc, 'scale': scale}
  824. >>> res = stats.goodness_of_fit(stats.norm, x, known_params=known_params,
  825. ... statistic='ks', random_state=rng)
  826. >>> res.statistic, res.pvalue
  827. (0.1119257570456813, 0.2788)
  828. The statistic matches exactly, but the p-value is estimated by forming
  829. a "Monte Carlo null distribution", that is, by explicitly drawing random
  830. samples from `scipy.stats.norm` with the provided parameters and
  831. calculating the stastic for each. The fraction of these statistic values
  832. at least as extreme as ``res.statistic`` approximates the exact p-value
  833. calculated by `scipy.stats.ks_1samp`.
  834. However, in many cases, we would prefer to test only that the data were
  835. sampled from one of *any* member of the normal distribution family, not
  836. specifically from the normal distribution with the location and scale
  837. fitted to the observed sample. In this case, Lilliefors [6]_ argued that
  838. the KS test is far too conservative (that is, the p-value overstates
  839. the actual probability of rejecting a true null hypothesis) and thus lacks
  840. power - the ability to reject the null hypothesis when the null hypothesis
  841. is actually false.
  842. Indeed, our p-value above is approximately 0.28, which is far too large
  843. to reject the null hypothesis at any common significance level.
  844. Consider why this might be. Note that in the KS test above, the statistic
  845. always compares data against the CDF of a normal distribution fitted to the
  846. *observed data*. This tends to reduce the value of the statistic for the
  847. observed data, but it is "unfair" when computing the statistic for other
  848. samples, such as those we randomly draw to form the Monte Carlo null
  849. distribution. It is easy to correct for this: whenever we compute the KS
  850. statistic of a sample, we use the CDF of a normal distribution fitted
  851. to *that sample*. The null distribution in this case has not been
  852. calculated exactly and is tyically approximated using Monte Carlo methods
  853. as described above. This is where `goodness_of_fit` excels.
  854. >>> res = stats.goodness_of_fit(stats.norm, x, statistic='ks',
  855. ... random_state=rng)
  856. >>> res.statistic, res.pvalue
  857. (0.1119257570456813, 0.0196)
  858. Indeed, this p-value is much smaller, and small enough to (correctly)
  859. reject the null hypothesis at common signficance levels, including 5% and
  860. 2.5%.
  861. However, the KS statistic is not very sensitive to all deviations from
  862. normality. The original advantage of the KS statistic was the ability
  863. to compute the null distribution theoretically, but a more sensitive
  864. statistic - resulting in a higher test power - can be used now that we can
  865. approximate the null distribution
  866. computationally. The Anderson-Darling statistic [1]_ tends to be more
  867. sensitive, and critical values of the this statistic have been tabulated
  868. for various significance levels and sample sizes using Monte Carlo methods.
  869. >>> res = stats.anderson(x, 'norm')
  870. >>> print(res.statistic)
  871. 1.2139573337497467
  872. >>> print(res.critical_values)
  873. [0.549 0.625 0.75 0.875 1.041]
  874. >>> print(res.significance_level)
  875. [15. 10. 5. 2.5 1. ]
  876. Here, the observed value of the statistic exceeds the critical value
  877. corresponding with a 1% significance level. This tells us that the p-value
  878. of the observed data is less than 1%, but what is it? We could interpolate
  879. from these (already-interpolated) values, but `goodness_of_fit` can
  880. estimate it directly.
  881. >>> res = stats.goodness_of_fit(stats.norm, x, statistic='ad',
  882. ... random_state=rng)
  883. >>> res.statistic, res.pvalue
  884. (1.2139573337497467, 0.0034)
  885. A further advantage is that use of `goodness_of_fit` is not limited to
  886. a particular set of distributions or conditions on which parameters
  887. are known versus which must be estimated from data. Instead,
  888. `goodness_of_fit` can estimate p-values relatively quickly for any
  889. distribution with a sufficiently fast and reliable ``fit`` method. For
  890. instance, here we perform a goodness of fit test using the Cramer-von Mises
  891. statistic against the Rayleigh distribution with known location and unknown
  892. scale.
  893. >>> rng = np.random.default_rng()
  894. >>> x = stats.chi(df=2.2, loc=0, scale=2).rvs(size=1000, random_state=rng)
  895. >>> res = stats.goodness_of_fit(stats.rayleigh, x, statistic='cvm',
  896. ... known_params={'loc': 0}, random_state=rng)
  897. This executes fairly quickly, but to check the reliability of the ``fit``
  898. method, we should inspect the fit result.
  899. >>> res.fit_result # location is as specified, and scale is reasonable
  900. params: FitParams(loc=0.0, scale=2.1026719844231243)
  901. success: True
  902. message: 'The fit was performed successfully.'
  903. >>> import matplotlib.pyplot as plt # matplotlib must be installed to plot
  904. >>> res.fit_result.plot()
  905. >>> plt.show()
  906. If the distribution is not fit to the observed data as well as possible,
  907. the test may not control the type I error rate, that is, the chance of
  908. rejecting the null hypothesis even when it is true.
  909. We should also look for extreme outliers in the null distribution that
  910. may be caused by unreliable fitting. These do not necessarily invalidate
  911. the result, but they tend to reduce the test's power.
  912. >>> _, ax = plt.subplots()
  913. >>> ax.hist(np.log10(res.null_distribution))
  914. >>> ax.set_xlabel("log10 of CVM statistic under the null hypothesis")
  915. >>> ax.set_ylabel("Frequency")
  916. >>> ax.set_title("Histogram of the Monte Carlo null distribution")
  917. >>> plt.show()
  918. This plot seems reassuring.
  919. If ``fit`` method is working reliably, and if the distribution of the test
  920. statistic is not particularly sensitive to the values of the fitted
  921. parameters, then the p-value provided by `goodness_of_fit` is expected to
  922. be a good approximation.
  923. >>> res.statistic, res.pvalue
  924. (0.2231991510248692, 0.0525)
  925. """
  926. args = _gof_iv(dist, data, known_params, fit_params, guessed_params,
  927. statistic, n_mc_samples, random_state)
  928. (dist, data, fixed_nhd_params, fixed_rfd_params, guessed_nhd_params,
  929. guessed_rfd_params, statistic, n_mc_samples_int, random_state) = args
  930. # Fit null hypothesis distribution to data
  931. nhd_fit_fun = _get_fit_fun(dist, data, guessed_nhd_params,
  932. fixed_nhd_params)
  933. nhd_vals = nhd_fit_fun(data)
  934. nhd_dist = dist(*nhd_vals)
  935. def rvs(size):
  936. return nhd_dist.rvs(size=size, random_state=random_state)
  937. # Define statistic
  938. fit_fun = _get_fit_fun(dist, data, guessed_rfd_params, fixed_rfd_params)
  939. compare_fun = _compare_dict[statistic]
  940. def statistic_fun(data, axis=-1):
  941. # Make things simple by always working along the last axis.
  942. data = np.moveaxis(data, axis, -1)
  943. rfd_vals = fit_fun(data)
  944. rfd_dist = dist(*rfd_vals)
  945. return compare_fun(rfd_dist, data)
  946. res = stats.monte_carlo_test(data, rvs, statistic_fun, vectorized=True,
  947. n_resamples=n_mc_samples, axis=-1,
  948. alternative='greater')
  949. opt_res = optimize.OptimizeResult()
  950. opt_res.success = True
  951. opt_res.message = "The fit was performed successfully."
  952. opt_res.x = nhd_vals
  953. # Only continuous distributions for now, hence discrete=False
  954. # There's no fundamental limitation; it's just that we're not using
  955. # stats.fit, discrete distributions don't have `fit` method, and
  956. # we haven't written any vectorized fit functions for a discrete
  957. # distribution yet.
  958. return GoodnessOfFitResult(FitResult(dist, data, False, opt_res),
  959. res.statistic, res.pvalue,
  960. res.null_distribution)
  961. def _get_fit_fun(dist, data, guessed_params, fixed_params):
  962. shape_names = [] if dist.shapes is None else dist.shapes.split(", ")
  963. param_names = shape_names + ['loc', 'scale']
  964. fparam_names = ['f'+name for name in param_names]
  965. all_fixed = not set(fparam_names).difference(fixed_params)
  966. guessed_shapes = [guessed_params.pop(x, None)
  967. for x in shape_names if x in guessed_params]
  968. # Define statistic, including fitting distribution to data
  969. if dist in _fit_funs:
  970. def fit_fun(data):
  971. params = _fit_funs[dist](data, **fixed_params)
  972. params = np.asarray(np.broadcast_arrays(*params))
  973. if params.ndim > 1:
  974. params = params[..., np.newaxis]
  975. return params
  976. elif all_fixed:
  977. def fit_fun(data):
  978. return [fixed_params[name] for name in fparam_names]
  979. else:
  980. def fit_fun_1d(data):
  981. return dist.fit(data, *guessed_shapes, **guessed_params,
  982. **fixed_params)
  983. def fit_fun(data):
  984. params = np.apply_along_axis(fit_fun_1d, axis=-1, arr=data)
  985. if params.ndim > 1:
  986. params = params.T[..., np.newaxis]
  987. return params
  988. return fit_fun
  989. # Vectorized fitting functions. These are to accept ND `data` in which each
  990. # row (slice along last axis) is a sample to fit and scalar fixed parameters.
  991. # They return a tuple of shape parameter arrays, each of shape data.shape[:-1].
  992. def _fit_norm(data, floc=None, fscale=None):
  993. loc = floc
  994. scale = fscale
  995. if loc is None and scale is None:
  996. loc = np.mean(data, axis=-1)
  997. scale = np.std(data, ddof=1, axis=-1)
  998. elif loc is None:
  999. loc = np.mean(data, axis=-1)
  1000. elif scale is None:
  1001. scale = np.sqrt(((data - loc)**2).mean(axis=-1))
  1002. return loc, scale
  1003. _fit_funs = {stats.norm: _fit_norm} # type: ignore[attr-defined]
  1004. # Vectorized goodness of fit statistic functions. These accept a frozen
  1005. # distribution object and `data` in which each row (slice along last axis) is
  1006. # a sample.
  1007. def _anderson_darling(dist, data):
  1008. x = np.sort(data, axis=-1)
  1009. n = data.shape[-1]
  1010. i = np.arange(1, n+1)
  1011. Si = (2*i - 1)/n * (dist.logcdf(x) + dist.logsf(x[..., ::-1]))
  1012. S = np.sum(Si, axis=-1)
  1013. return -n - S
  1014. def _compute_dplus(cdfvals): # adapted from _stats_py before gh-17062
  1015. n = cdfvals.shape[-1]
  1016. return (np.arange(1.0, n + 1) / n - cdfvals).max(axis=-1)
  1017. def _compute_dminus(cdfvals, axis=-1):
  1018. n = cdfvals.shape[-1]
  1019. return (cdfvals - np.arange(0.0, n)/n).max(axis=-1)
  1020. def _kolmogorov_smirnov(dist, data):
  1021. x = np.sort(data, axis=-1)
  1022. cdfvals = dist.cdf(x)
  1023. Dplus = _compute_dplus(cdfvals) # always works along last axis
  1024. Dminus = _compute_dminus(cdfvals)
  1025. return np.maximum(Dplus, Dminus)
  1026. def _cramer_von_mises(dist, data):
  1027. x = np.sort(data, axis=-1)
  1028. n = data.shape[-1]
  1029. cdfvals = dist.cdf(x)
  1030. u = (2*np.arange(1, n+1) - 1)/(2*n)
  1031. w = 1 / (12*n) + np.sum((u - cdfvals)**2, axis=-1)
  1032. return w
  1033. _compare_dict = {"ad": _anderson_darling, "ks": _kolmogorov_smirnov,
  1034. "cvm": _cramer_von_mises}
  1035. def _gof_iv(dist, data, known_params, fit_params, guessed_params, statistic,
  1036. n_mc_samples, random_state):
  1037. if not isinstance(dist, stats.rv_continuous):
  1038. message = ("`dist` must be a (non-frozen) instance of "
  1039. "`stats.rv_continuous`.")
  1040. raise TypeError(message)
  1041. data = np.asarray(data, dtype=float)
  1042. if not data.ndim == 1:
  1043. message = "`data` must be a one-dimensional array of numbers."
  1044. raise ValueError(message)
  1045. # Leave validation of these key/value pairs to the `fit` method,
  1046. # but collect these into dictionaries that will be used
  1047. known_params = known_params or dict()
  1048. fit_params = fit_params or dict()
  1049. guessed_params = guessed_params or dict()
  1050. known_params_f = {("f"+key): val for key, val in known_params.items()}
  1051. fit_params_f = {("f"+key): val for key, val in fit_params.items()}
  1052. # These the the values of parameters of the null distribution family
  1053. # with which resamples are drawn
  1054. fixed_nhd_params = known_params_f.copy()
  1055. fixed_nhd_params.update(fit_params_f)
  1056. # These are fixed when fitting the distribution family to resamples
  1057. fixed_rfd_params = known_params_f.copy()
  1058. # These are used as guesses when fitting the distribution family to
  1059. # the original data
  1060. guessed_nhd_params = guessed_params.copy()
  1061. # These are used as guesses when fitting the distribution family to
  1062. # resamples
  1063. guessed_rfd_params = fit_params.copy()
  1064. guessed_rfd_params.update(guessed_params)
  1065. statistics = {'ad', 'ks', 'cvm'}
  1066. if statistic.lower() not in statistics:
  1067. message = f"`statistic` must be one of {statistics}."
  1068. raise ValueError(message)
  1069. n_mc_samples_int = int(n_mc_samples)
  1070. if n_mc_samples_int != n_mc_samples:
  1071. message = "`n_mc_samples` must be an integer."
  1072. raise TypeError(message)
  1073. random_state = check_random_state(random_state)
  1074. return (dist, data, fixed_nhd_params, fixed_rfd_params, guessed_nhd_params,
  1075. guessed_rfd_params, statistic, n_mc_samples_int, random_state)