_qmc.py 105 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735173617371738173917401741174217431744174517461747174817491750175117521753175417551756175717581759176017611762176317641765176617671768176917701771177217731774177517761777177817791780178117821783178417851786178717881789179017911792179317941795179617971798179918001801180218031804180518061807180818091810181118121813181418151816181718181819182018211822182318241825182618271828182918301831183218331834183518361837183818391840184118421843184418451846184718481849185018511852185318541855185618571858185918601861186218631864186518661867186818691870187118721873187418751876187718781879188018811882188318841885188618871888188918901891189218931894189518961897189818991900190119021903190419051906190719081909191019111912191319141915191619171918191919201921192219231924192519261927192819291930193119321933193419351936193719381939194019411942194319441945194619471948194919501951195219531954195519561957195819591960196119621963196419651966196719681969197019711972197319741975197619771978197919801981198219831984198519861987198819891990199119921993199419951996199719981999200020012002200320042005200620072008200920102011201220132014201520162017201820192020202120222023202420252026202720282029203020312032203320342035203620372038203920402041204220432044204520462047204820492050205120522053205420552056205720582059206020612062206320642065206620672068206920702071207220732074207520762077207820792080208120822083208420852086208720882089209020912092209320942095209620972098209921002101210221032104210521062107210821092110211121122113211421152116211721182119212021212122212321242125212621272128212921302131213221332134213521362137213821392140214121422143214421452146214721482149215021512152215321542155215621572158215921602161216221632164216521662167216821692170217121722173217421752176217721782179218021812182218321842185218621872188218921902191219221932194219521962197219821992200220122022203220422052206220722082209221022112212221322142215221622172218221922202221222222232224222522262227222822292230223122322233223422352236223722382239224022412242224322442245224622472248224922502251225222532254225522562257225822592260226122622263226422652266226722682269227022712272227322742275227622772278227922802281228222832284228522862287228822892290229122922293229422952296229722982299230023012302230323042305230623072308230923102311231223132314231523162317231823192320232123222323232423252326232723282329233023312332233323342335233623372338233923402341234223432344234523462347234823492350235123522353235423552356235723582359236023612362236323642365236623672368236923702371237223732374237523762377237823792380238123822383238423852386238723882389239023912392239323942395239623972398239924002401240224032404240524062407240824092410241124122413241424152416241724182419242024212422242324242425242624272428242924302431243224332434243524362437243824392440244124422443244424452446244724482449245024512452245324542455245624572458245924602461246224632464246524662467246824692470247124722473247424752476247724782479248024812482248324842485248624872488248924902491249224932494249524962497249824992500250125022503250425052506250725082509251025112512251325142515251625172518251925202521252225232524252525262527252825292530253125322533253425352536253725382539254025412542254325442545254625472548254925502551255225532554255525562557255825592560256125622563256425652566256725682569257025712572257325742575257625772578257925802581258225832584258525862587258825892590259125922593259425952596259725982599260026012602260326042605260626072608260926102611261226132614261526162617261826192620262126222623262426252626262726282629263026312632263326342635263626372638263926402641264226432644264526462647264826492650265126522653265426552656265726582659266026612662266326642665266626672668266926702671267226732674267526762677267826792680268126822683268426852686268726882689269026912692269326942695269626972698269927002701270227032704270527062707270827092710271127122713271427152716271727182719272027212722272327242725272627272728272927302731273227332734273527362737273827392740274127422743274427452746274727482749275027512752275327542755275627572758275927602761276227632764276527662767276827692770277127722773277427752776277727782779278027812782278327842785278627872788278927902791279227932794279527962797279827992800280128022803280428052806280728082809281028112812281328142815281628172818281928202821282228232824282528262827282828292830283128322833283428352836283728382839284028412842284328442845284628472848284928502851285228532854285528562857285828592860286128622863286428652866286728682869287028712872287328742875287628772878287928802881288228832884288528862887288828892890289128922893289428952896289728982899290029012902290329042905290629072908290929102911291229132914291529162917291829192920292129222923292429252926292729282929293029312932293329342935293629372938293929402941294229432944294529462947294829492950295129522953295429552956
  1. """Quasi-Monte Carlo engines and helpers."""
  2. import copy
  3. import math
  4. import numbers
  5. import os
  6. import warnings
  7. from abc import ABC, abstractmethod
  8. from functools import partial
  9. from typing import (
  10. ClassVar,
  11. Literal,
  12. overload,
  13. TYPE_CHECKING,
  14. )
  15. from collections.abc import Callable
  16. import numpy as np
  17. from scipy._lib._util import DecimalNumber, GeneratorType, IntNumber, SeedType
  18. if TYPE_CHECKING:
  19. import numpy.typing as npt
  20. import scipy.stats as stats
  21. from scipy._lib._util import rng_integers, _rng_spawn, _transition_to_rng
  22. from scipy.sparse.csgraph import minimum_spanning_tree
  23. from scipy.spatial import distance, Voronoi
  24. from scipy.special import gammainc
  25. from ._sobol import (
  26. _initialize_v, _cscramble, _fill_p_cumulative, _draw, _fast_forward,
  27. _categorize, _MAXDIM
  28. )
  29. from ._qmc_cy import (
  30. _cy_wrapper_centered_discrepancy,
  31. _cy_wrapper_wrap_around_discrepancy,
  32. _cy_wrapper_mixture_discrepancy,
  33. _cy_wrapper_l2_star_discrepancy,
  34. _cy_wrapper_update_discrepancy,
  35. _cy_van_der_corput_scrambled,
  36. _cy_van_der_corput,
  37. )
  38. __all__ = ['scale', 'discrepancy', 'geometric_discrepancy', 'update_discrepancy',
  39. 'QMCEngine', 'Sobol', 'Halton', 'LatinHypercube', 'PoissonDisk',
  40. 'MultinomialQMC', 'MultivariateNormalQMC']
  41. @overload
  42. def check_random_state(seed: IntNumber | None = ...) -> np.random.Generator:
  43. ...
  44. @overload
  45. def check_random_state(seed: GeneratorType) -> GeneratorType:
  46. ...
  47. # Based on scipy._lib._util.check_random_state
  48. # This is going to be removed at the end of the SPEC 7 transition,
  49. # so I'll just leave the argument name `seed` alone
  50. def check_random_state(seed=None):
  51. """Turn `seed` into a `numpy.random.Generator` instance.
  52. Parameters
  53. ----------
  54. seed : {None, int, `numpy.random.Generator`, `numpy.random.RandomState`}, optional
  55. If `seed` is an int or None, a new `numpy.random.Generator` is
  56. created using ``np.random.default_rng(seed)``.
  57. If `seed` is already a ``Generator`` or ``RandomState`` instance, then
  58. the provided instance is used.
  59. Returns
  60. -------
  61. seed : {`numpy.random.Generator`, `numpy.random.RandomState`}
  62. Random number generator.
  63. """
  64. if seed is None or isinstance(seed, numbers.Integral | np.integer):
  65. return np.random.default_rng(seed)
  66. elif isinstance(seed, np.random.RandomState | np.random.Generator):
  67. return seed
  68. else:
  69. raise ValueError(f'{seed!r} cannot be used to seed a'
  70. ' numpy.random.Generator instance')
  71. def scale(
  72. sample: "npt.ArrayLike",
  73. l_bounds: "npt.ArrayLike",
  74. u_bounds: "npt.ArrayLike",
  75. *,
  76. reverse: bool = False
  77. ) -> np.ndarray:
  78. r"""Sample scaling from unit hypercube to different bounds.
  79. To convert a sample from :math:`[0, 1)` to :math:`[a, b), b>a`,
  80. with :math:`a` the lower bounds and :math:`b` the upper bounds.
  81. The following transformation is used:
  82. .. math::
  83. (b - a) \cdot \text{sample} + a
  84. Parameters
  85. ----------
  86. sample : array_like (n, d)
  87. Sample to scale.
  88. l_bounds, u_bounds : array_like (d,)
  89. Lower and upper bounds (resp. :math:`a`, :math:`b`) of transformed
  90. data. If `reverse` is True, range of the original data to transform
  91. to the unit hypercube.
  92. reverse : bool, optional
  93. Reverse the transformation from different bounds to the unit hypercube.
  94. Default is False.
  95. Returns
  96. -------
  97. sample : array_like (n, d)
  98. Scaled sample.
  99. Examples
  100. --------
  101. Transform 3 samples in the unit hypercube to bounds:
  102. >>> from scipy.stats import qmc
  103. >>> l_bounds = [-2, 0]
  104. >>> u_bounds = [6, 5]
  105. >>> sample = [[0.5 , 0.75],
  106. ... [0.5 , 0.5],
  107. ... [0.75, 0.25]]
  108. >>> sample_scaled = qmc.scale(sample, l_bounds, u_bounds)
  109. >>> sample_scaled
  110. array([[2. , 3.75],
  111. [2. , 2.5 ],
  112. [4. , 1.25]])
  113. And convert back to the unit hypercube:
  114. >>> sample_ = qmc.scale(sample_scaled, l_bounds, u_bounds, reverse=True)
  115. >>> sample_
  116. array([[0.5 , 0.75],
  117. [0.5 , 0.5 ],
  118. [0.75, 0.25]])
  119. """
  120. sample = np.asarray(sample)
  121. # Checking bounds and sample
  122. if not sample.ndim == 2:
  123. raise ValueError('Sample is not a 2D array')
  124. lower, upper = _validate_bounds(
  125. l_bounds=l_bounds, u_bounds=u_bounds, d=sample.shape[1]
  126. )
  127. if not reverse:
  128. # Checking that sample is within the hypercube
  129. if (sample.max() > 1.) or (sample.min() < 0.):
  130. raise ValueError('Sample is not in unit hypercube')
  131. return sample * (upper - lower) + lower
  132. else:
  133. # Checking that sample is within the bounds
  134. if not (np.all(sample >= lower) and np.all(sample <= upper)):
  135. raise ValueError('Sample is out of bounds')
  136. return (sample - lower) / (upper - lower)
  137. def _ensure_in_unit_hypercube(sample: "npt.ArrayLike") -> np.ndarray:
  138. """Ensure that sample is a 2D array and is within a unit hypercube
  139. Parameters
  140. ----------
  141. sample : array_like (n, d)
  142. A 2D array of points.
  143. Returns
  144. -------
  145. np.ndarray
  146. The array interpretation of the input sample
  147. Raises
  148. ------
  149. ValueError
  150. If the input is not a 2D array or contains points outside of
  151. a unit hypercube.
  152. """
  153. sample = np.asarray(sample, dtype=np.float64, order="C")
  154. if not sample.ndim == 2:
  155. raise ValueError("Sample is not a 2D array")
  156. if (sample.max() > 1.) or (sample.min() < 0.):
  157. raise ValueError("Sample is not in unit hypercube")
  158. return sample
  159. def discrepancy(
  160. sample: "npt.ArrayLike",
  161. *,
  162. iterative: bool = False,
  163. method: Literal["CD", "WD", "MD", "L2-star"] = "CD",
  164. workers: IntNumber = 1) -> float:
  165. """Discrepancy of a given sample.
  166. Parameters
  167. ----------
  168. sample : array_like (n, d)
  169. The sample to compute the discrepancy from.
  170. iterative : bool, optional
  171. Must be False if not using it for updating the discrepancy.
  172. Default is False. Refer to the notes for more details.
  173. method : str, optional
  174. Type of discrepancy, can be ``CD``, ``WD``, ``MD`` or ``L2-star``.
  175. Refer to the notes for more details. Default is ``CD``.
  176. workers : int, optional
  177. Number of workers to use for parallel processing. If -1 is given all
  178. CPU threads are used. Default is 1.
  179. Returns
  180. -------
  181. discrepancy : float
  182. Discrepancy.
  183. See Also
  184. --------
  185. geometric_discrepancy
  186. Notes
  187. -----
  188. The discrepancy is a uniformity criterion used to assess the space filling
  189. of a number of samples in a hypercube. A discrepancy quantifies the
  190. distance between the continuous uniform distribution on a hypercube and the
  191. discrete uniform distribution on :math:`n` distinct sample points.
  192. The lower the value is, the better the coverage of the parameter space is.
  193. For a collection of subsets of the hypercube, the discrepancy is the
  194. difference between the fraction of sample points in one of those
  195. subsets and the volume of that subset. There are different definitions of
  196. discrepancy corresponding to different collections of subsets. Some
  197. versions take a root mean square difference over subsets instead of
  198. a maximum.
  199. A measure of uniformity is reasonable if it satisfies the following
  200. criteria [1]_:
  201. 1. It is invariant under permuting factors and/or runs.
  202. 2. It is invariant under rotation of the coordinates.
  203. 3. It can measure not only uniformity of the sample over the hypercube,
  204. but also the projection uniformity of the sample over non-empty
  205. subset of lower dimension hypercubes.
  206. 4. There is some reasonable geometric meaning.
  207. 5. It is easy to compute.
  208. 6. It satisfies the Koksma-Hlawka-like inequality.
  209. 7. It is consistent with other criteria in experimental design.
  210. Four methods are available:
  211. * ``CD``: Centered Discrepancy - subspace involves a corner of the
  212. hypercube
  213. * ``WD``: Wrap-around Discrepancy - subspace can wrap around bounds
  214. * ``MD``: Mixture Discrepancy - mix between CD/WD covering more criteria
  215. * ``L2-star``: L2-star discrepancy - like CD BUT variant to rotation
  216. Methods ``CD``, ``WD``, and ``MD`` implement the right hand side of equations
  217. 9, 10, and 18 of [2]_, respectively; the square root is not taken. On the
  218. other hand, ``L2-star`` computes the quantity given by equation 10 of
  219. [3]_ as implemented by subsequent equations; the square root is taken.
  220. Lastly, using ``iterative=True``, it is possible to compute the
  221. discrepancy as if we had :math:`n+1` samples. This is useful if we want
  222. to add a point to a sampling and check the candidate which would give the
  223. lowest discrepancy. Then you could just update the discrepancy with
  224. each candidate using `update_discrepancy`. This method is faster than
  225. computing the discrepancy for a large number of candidates.
  226. References
  227. ----------
  228. .. [1] Fang et al. "Design and modeling for computer experiments".
  229. Computer Science and Data Analysis Series, 2006.
  230. .. [2] Zhou Y.-D. et al. "Mixture discrepancy for quasi-random point sets."
  231. Journal of Complexity, 29 (3-4) , pp. 283-301, 2013.
  232. .. [3] T. T. Warnock. "Computational investigations of low discrepancy
  233. point sets." Applications of Number Theory to Numerical
  234. Analysis, Academic Press, pp. 319-343, 1972.
  235. Examples
  236. --------
  237. Calculate the quality of the sample using the discrepancy:
  238. >>> import numpy as np
  239. >>> from scipy.stats import qmc
  240. >>> space = np.array([[1, 3], [2, 6], [3, 2], [4, 5], [5, 1], [6, 4]])
  241. >>> l_bounds = [0.5, 0.5]
  242. >>> u_bounds = [6.5, 6.5]
  243. >>> space = qmc.scale(space, l_bounds, u_bounds, reverse=True)
  244. >>> space
  245. array([[0.08333333, 0.41666667],
  246. [0.25 , 0.91666667],
  247. [0.41666667, 0.25 ],
  248. [0.58333333, 0.75 ],
  249. [0.75 , 0.08333333],
  250. [0.91666667, 0.58333333]])
  251. >>> qmc.discrepancy(space)
  252. 0.008142039609053464
  253. We can also compute iteratively the ``CD`` discrepancy by using
  254. ``iterative=True``.
  255. >>> disc_init = qmc.discrepancy(space[:-1], iterative=True)
  256. >>> disc_init
  257. 0.04769081147119336
  258. >>> qmc.update_discrepancy(space[-1], space[:-1], disc_init)
  259. 0.008142039609053513
  260. """
  261. sample = _ensure_in_unit_hypercube(sample)
  262. workers = _validate_workers(workers)
  263. methods = {
  264. "CD": _cy_wrapper_centered_discrepancy,
  265. "WD": _cy_wrapper_wrap_around_discrepancy,
  266. "MD": _cy_wrapper_mixture_discrepancy,
  267. "L2-star": _cy_wrapper_l2_star_discrepancy,
  268. }
  269. if method in methods:
  270. return methods[method](sample, iterative, workers=workers)
  271. else:
  272. raise ValueError(f"{method!r} is not a valid method. It must be one of"
  273. f" {set(methods)!r}")
  274. def geometric_discrepancy(
  275. sample: "npt.ArrayLike",
  276. method: Literal["mindist", "mst"] = "mindist",
  277. metric: str = "euclidean") -> float:
  278. """Discrepancy of a given sample based on its geometric properties.
  279. Parameters
  280. ----------
  281. sample : array_like (n, d)
  282. The sample to compute the discrepancy from.
  283. method : {"mindist", "mst"}, optional
  284. The method to use. One of ``mindist`` for minimum distance (default)
  285. or ``mst`` for minimum spanning tree.
  286. metric : str or callable, optional
  287. The distance metric to use. See the documentation
  288. for `scipy.spatial.distance.pdist` for the available metrics and
  289. the default.
  290. Returns
  291. -------
  292. discrepancy : float
  293. Discrepancy (higher values correspond to greater sample uniformity).
  294. See Also
  295. --------
  296. discrepancy
  297. Notes
  298. -----
  299. The discrepancy can serve as a simple measure of quality of a random sample.
  300. This measure is based on the geometric properties of the distribution of points
  301. in the sample, such as the minimum distance between any pair of points, or
  302. the mean edge length in a minimum spanning tree.
  303. The higher the value is, the better the coverage of the parameter space is.
  304. Note that this is different from `scipy.stats.qmc.discrepancy`, where lower
  305. values correspond to higher quality of the sample.
  306. Also note that when comparing different sampling strategies using this function,
  307. the sample size must be kept constant.
  308. It is possible to calculate two metrics from the minimum spanning tree:
  309. the mean edge length and the standard deviation of edges lengths. Using
  310. both metrics offers a better picture of uniformity than either metric alone,
  311. with higher mean and lower standard deviation being preferable (see [1]_
  312. for a brief discussion). This function currently only calculates the mean
  313. edge length.
  314. References
  315. ----------
  316. .. [1] Franco J. et al. "Minimum Spanning Tree: A new approach to assess the quality
  317. of the design of computer experiments." Chemometrics and Intelligent Laboratory
  318. Systems, 97 (2), pp. 164-169, 2009.
  319. Examples
  320. --------
  321. Calculate the quality of the sample using the minimum euclidean distance
  322. (the defaults):
  323. >>> import numpy as np
  324. >>> from scipy.stats import qmc
  325. >>> rng = np.random.default_rng(191468432622931918890291693003068437394)
  326. >>> sample = qmc.LatinHypercube(d=2, rng=rng).random(50)
  327. >>> qmc.geometric_discrepancy(sample)
  328. 0.03708161435687876
  329. Calculate the quality using the mean edge length in the minimum
  330. spanning tree:
  331. >>> qmc.geometric_discrepancy(sample, method='mst')
  332. 0.1105149978798376
  333. Display the minimum spanning tree and the points with
  334. the smallest distance:
  335. >>> import matplotlib.pyplot as plt
  336. >>> from matplotlib.lines import Line2D
  337. >>> from scipy.sparse.csgraph import minimum_spanning_tree
  338. >>> from scipy.spatial.distance import pdist, squareform
  339. >>> dist = pdist(sample)
  340. >>> mst = minimum_spanning_tree(squareform(dist))
  341. >>> edges = np.where(mst.toarray() > 0)
  342. >>> edges = np.asarray(edges).T
  343. >>> min_dist = np.min(dist)
  344. >>> min_idx = np.argwhere(squareform(dist) == min_dist)[0]
  345. >>> fig, ax = plt.subplots(figsize=(10, 5))
  346. >>> _ = ax.set(aspect='equal', xlabel=r'$x_1$', ylabel=r'$x_2$',
  347. ... xlim=[0, 1], ylim=[0, 1])
  348. >>> for edge in edges:
  349. ... ax.plot(sample[edge, 0], sample[edge, 1], c='k')
  350. >>> ax.scatter(sample[:, 0], sample[:, 1])
  351. >>> ax.add_patch(plt.Circle(sample[min_idx[0]], min_dist, color='red', fill=False))
  352. >>> markers = [
  353. ... Line2D([0], [0], marker='o', lw=0, label='Sample points'),
  354. ... Line2D([0], [0], color='k', label='Minimum spanning tree'),
  355. ... Line2D([0], [0], marker='o', lw=0, markerfacecolor='w', markeredgecolor='r',
  356. ... label='Minimum point-to-point distance'),
  357. ... ]
  358. >>> ax.legend(handles=markers, loc='center left', bbox_to_anchor=(1, 0.5));
  359. >>> plt.show()
  360. """
  361. sample = _ensure_in_unit_hypercube(sample)
  362. if sample.shape[0] < 2:
  363. raise ValueError("Sample must contain at least two points")
  364. distances = distance.pdist(sample, metric=metric) # type: ignore[call-overload]
  365. if np.any(distances == 0.0):
  366. warnings.warn("Sample contains duplicate points.", stacklevel=2)
  367. if method == "mindist":
  368. return np.min(distances[distances.nonzero()])
  369. elif method == "mst":
  370. fully_connected_graph = distance.squareform(distances)
  371. mst = minimum_spanning_tree(fully_connected_graph)
  372. distances = mst[mst.nonzero()]
  373. # TODO consider returning both the mean and the standard deviation
  374. # see [1] for a discussion
  375. return np.mean(distances)
  376. else:
  377. raise ValueError(f"{method!r} is not a valid method. "
  378. f"It must be one of {{'mindist', 'mst'}}")
  379. def update_discrepancy(
  380. x_new: "npt.ArrayLike",
  381. sample: "npt.ArrayLike",
  382. initial_disc: DecimalNumber) -> float:
  383. """Update the centered discrepancy with a new sample.
  384. Parameters
  385. ----------
  386. x_new : array_like (1, d)
  387. The new sample to add in `sample`.
  388. sample : array_like (n, d)
  389. The initial sample.
  390. initial_disc : float
  391. Centered discrepancy of the `sample`.
  392. Returns
  393. -------
  394. discrepancy : float
  395. Centered discrepancy of the sample composed of `x_new` and `sample`.
  396. Examples
  397. --------
  398. We can also compute iteratively the discrepancy by using
  399. ``iterative=True``.
  400. >>> import numpy as np
  401. >>> from scipy.stats import qmc
  402. >>> space = np.array([[1, 3], [2, 6], [3, 2], [4, 5], [5, 1], [6, 4]])
  403. >>> l_bounds = [0.5, 0.5]
  404. >>> u_bounds = [6.5, 6.5]
  405. >>> space = qmc.scale(space, l_bounds, u_bounds, reverse=True)
  406. >>> disc_init = qmc.discrepancy(space[:-1], iterative=True)
  407. >>> disc_init
  408. 0.04769081147119336
  409. >>> qmc.update_discrepancy(space[-1], space[:-1], disc_init)
  410. 0.008142039609053513
  411. """
  412. sample = np.asarray(sample, dtype=np.float64, order="C")
  413. x_new = np.asarray(x_new, dtype=np.float64, order="C")
  414. # Checking that sample is within the hypercube and 2D
  415. if not sample.ndim == 2:
  416. raise ValueError('Sample is not a 2D array')
  417. if (sample.max() > 1.) or (sample.min() < 0.):
  418. raise ValueError('Sample is not in unit hypercube')
  419. # Checking that x_new is within the hypercube and 1D
  420. if not x_new.ndim == 1:
  421. raise ValueError('x_new is not a 1D array')
  422. if not (np.all(x_new >= 0) and np.all(x_new <= 1)):
  423. raise ValueError('x_new is not in unit hypercube')
  424. if x_new.shape[0] != sample.shape[1]:
  425. raise ValueError("x_new and sample must be broadcastable")
  426. return _cy_wrapper_update_discrepancy(x_new, sample, initial_disc)
  427. def _perturb_discrepancy(sample: np.ndarray, i1: int, i2: int, k: int,
  428. disc: float):
  429. """Centered discrepancy after an elementary perturbation of a LHS.
  430. An elementary perturbation consists of an exchange of coordinates between
  431. two points: ``sample[i1, k] <-> sample[i2, k]``. By construction,
  432. this operation conserves the LHS properties.
  433. Parameters
  434. ----------
  435. sample : array_like (n, d)
  436. The sample (before permutation) to compute the discrepancy from.
  437. i1 : int
  438. The first line of the elementary permutation.
  439. i2 : int
  440. The second line of the elementary permutation.
  441. k : int
  442. The column of the elementary permutation.
  443. disc : float
  444. Centered discrepancy of the design before permutation.
  445. Returns
  446. -------
  447. discrepancy : float
  448. Centered discrepancy of the design after permutation.
  449. References
  450. ----------
  451. .. [1] Jin et al. "An efficient algorithm for constructing optimal design
  452. of computer experiments", Journal of Statistical Planning and
  453. Inference, 2005.
  454. """
  455. n = sample.shape[0]
  456. z_ij = sample - 0.5
  457. # Eq (19)
  458. c_i1j = (1. / n ** 2.
  459. * np.prod(0.5 * (2. + abs(z_ij[i1, :])
  460. + abs(z_ij) - abs(z_ij[i1, :] - z_ij)), axis=1))
  461. c_i2j = (1. / n ** 2.
  462. * np.prod(0.5 * (2. + abs(z_ij[i2, :])
  463. + abs(z_ij) - abs(z_ij[i2, :] - z_ij)), axis=1))
  464. # Eq (20)
  465. c_i1i1 = (1. / n ** 2 * np.prod(1 + abs(z_ij[i1, :]))
  466. - 2. / n * np.prod(1. + 0.5 * abs(z_ij[i1, :])
  467. - 0.5 * z_ij[i1, :] ** 2))
  468. c_i2i2 = (1. / n ** 2 * np.prod(1 + abs(z_ij[i2, :]))
  469. - 2. / n * np.prod(1. + 0.5 * abs(z_ij[i2, :])
  470. - 0.5 * z_ij[i2, :] ** 2))
  471. # Eq (22), typo in the article in the denominator i2 -> i1
  472. num = (2 + abs(z_ij[i2, k]) + abs(z_ij[:, k])
  473. - abs(z_ij[i2, k] - z_ij[:, k]))
  474. denum = (2 + abs(z_ij[i1, k]) + abs(z_ij[:, k])
  475. - abs(z_ij[i1, k] - z_ij[:, k]))
  476. gamma = num / denum
  477. # Eq (23)
  478. c_p_i1j = gamma * c_i1j
  479. # Eq (24)
  480. c_p_i2j = c_i2j / gamma
  481. alpha = (1 + abs(z_ij[i2, k])) / (1 + abs(z_ij[i1, k]))
  482. beta = (2 - abs(z_ij[i2, k])) / (2 - abs(z_ij[i1, k]))
  483. g_i1 = np.prod(1. + abs(z_ij[i1, :]))
  484. g_i2 = np.prod(1. + abs(z_ij[i2, :]))
  485. h_i1 = np.prod(1. + 0.5 * abs(z_ij[i1, :]) - 0.5 * (z_ij[i1, :] ** 2))
  486. h_i2 = np.prod(1. + 0.5 * abs(z_ij[i2, :]) - 0.5 * (z_ij[i2, :] ** 2))
  487. # Eq (25), typo in the article g is missing
  488. c_p_i1i1 = ((g_i1 * alpha) / (n ** 2) - 2. * alpha * beta * h_i1 / n)
  489. # Eq (26), typo in the article n ** 2
  490. c_p_i2i2 = ((g_i2 / ((n ** 2) * alpha)) - (2. * h_i2 / (n * alpha * beta)))
  491. # Eq (26)
  492. sum_ = c_p_i1j - c_i1j + c_p_i2j - c_i2j
  493. mask = np.ones(n, dtype=bool)
  494. mask[[i1, i2]] = False
  495. sum_ = sum(sum_[mask])
  496. disc_ep = (disc + c_p_i1i1 - c_i1i1 + c_p_i2i2 - c_i2i2 + 2 * sum_)
  497. return disc_ep
  498. def primes_from_2_to(n: int) -> np.ndarray:
  499. """Prime numbers from 2 to *n*.
  500. Parameters
  501. ----------
  502. n : int
  503. Sup bound with ``n >= 6``.
  504. Returns
  505. -------
  506. primes : list(int)
  507. Primes in ``2 <= p < n``.
  508. Notes
  509. -----
  510. Taken from [1]_ by P.T. Roy, written consent given on 23.04.2021
  511. by the original author, Bruno Astrolino, for free use in SciPy under
  512. the 3-clause BSD.
  513. References
  514. ----------
  515. .. [1] `StackOverflow <https://stackoverflow.com/questions/2068372>`_.
  516. """
  517. sieve = np.ones(n // 3 + (n % 6 == 2), dtype=bool)
  518. for i in range(1, int(n ** 0.5) // 3 + 1):
  519. k = 3 * i + 1 | 1
  520. sieve[k * k // 3::2 * k] = False
  521. sieve[k * (k - 2 * (i & 1) + 4) // 3::2 * k] = False
  522. return np.r_[2, 3, ((3 * np.nonzero(sieve)[0][1:] + 1) | 1)]
  523. def n_primes(n: IntNumber) -> list[int]:
  524. """List of the n-first prime numbers.
  525. Parameters
  526. ----------
  527. n : int
  528. Number of prime numbers wanted.
  529. Returns
  530. -------
  531. primes : list(int)
  532. List of primes.
  533. """
  534. primes = [2, 3, 5, 7, 11, 13, 17, 19, 23, 29, 31, 37, 41, 43, 47, 53, 59,
  535. 61, 67, 71, 73, 79, 83, 89, 97, 101, 103, 107, 109, 113, 127,
  536. 131, 137, 139, 149, 151, 157, 163, 167, 173, 179, 181, 191, 193,
  537. 197, 199, 211, 223, 227, 229, 233, 239, 241, 251, 257, 263, 269,
  538. 271, 277, 281, 283, 293, 307, 311, 313, 317, 331, 337, 347, 349,
  539. 353, 359, 367, 373, 379, 383, 389, 397, 401, 409, 419, 421, 431,
  540. 433, 439, 443, 449, 457, 461, 463, 467, 479, 487, 491, 499, 503,
  541. 509, 521, 523, 541, 547, 557, 563, 569, 571, 577, 587, 593, 599,
  542. 601, 607, 613, 617, 619, 631, 641, 643, 647, 653, 659, 661, 673,
  543. 677, 683, 691, 701, 709, 719, 727, 733, 739, 743, 751, 757, 761,
  544. 769, 773, 787, 797, 809, 811, 821, 823, 827, 829, 839, 853, 857,
  545. 859, 863, 877, 881, 883, 887, 907, 911, 919, 929, 937, 941, 947,
  546. 953, 967, 971, 977, 983, 991, 997][:n]
  547. if len(primes) < n:
  548. big_number = 2000
  549. while 'Not enough primes':
  550. primes = primes_from_2_to(big_number)[:n] # type: ignore
  551. if len(primes) == n:
  552. break
  553. big_number += 1000
  554. return primes
  555. def _van_der_corput_permutations(
  556. base: IntNumber, *, rng: SeedType = None
  557. ) -> np.ndarray:
  558. """Permutations for scrambling a Van der Corput sequence.
  559. Parameters
  560. ----------
  561. base : int
  562. Base of the sequence.
  563. rng : `numpy.random.Generator`, optional
  564. Pseudorandom number generator state. When `rng` is None, a new
  565. `numpy.random.Generator` is created using entropy from the
  566. operating system. Types other than `numpy.random.Generator` are
  567. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  568. .. versionchanged:: 1.15.0
  569. As part of the `SPEC-007 <https://scientific-python.org/specs/spec-0007/>`_
  570. transition from use of `numpy.random.RandomState` to
  571. `numpy.random.Generator`, this keyword was changed from `seed` to
  572. `rng`. During the transition, the behavior documented above is not
  573. accurate; see `check_random_state` for actual behavior. After the
  574. transition, this admonition can be removed.
  575. Returns
  576. -------
  577. permutations : array_like
  578. Permutation indices.
  579. Notes
  580. -----
  581. In Algorithm 1 of Owen 2017, a permutation of `np.arange(base)` is
  582. created for each positive integer `k` such that ``1 - base**-k < 1``
  583. using floating-point arithmetic. For double precision floats, the
  584. condition ``1 - base**-k < 1`` can also be written as ``base**-k >
  585. 2**-54``, which makes it more apparent how many permutations we need
  586. to create.
  587. """
  588. rng = check_random_state(rng)
  589. count = math.ceil(54 / math.log2(base)) - 1
  590. permutations = np.repeat(np.arange(base)[None], count, axis=0)
  591. for perm in permutations:
  592. rng.shuffle(perm)
  593. return permutations
  594. def van_der_corput(
  595. n: IntNumber,
  596. base: IntNumber = 2,
  597. *,
  598. start_index: IntNumber = 0,
  599. scramble: bool = False,
  600. permutations: "npt.ArrayLike | None" = None,
  601. rng: SeedType = None,
  602. workers: IntNumber = 1) -> np.ndarray:
  603. """Van der Corput sequence.
  604. Pseudo-random number generator based on a b-adic expansion.
  605. Scrambling uses permutations of the remainders (see [1]_). Multiple
  606. permutations are applied to construct a point. The sequence of
  607. permutations has to be the same for all points of the sequence.
  608. Parameters
  609. ----------
  610. n : int
  611. Number of element of the sequence.
  612. base : int, optional
  613. Base of the sequence. Default is 2.
  614. start_index : int, optional
  615. Index to start the sequence from. Default is 0.
  616. scramble : bool, optional
  617. If True, use Owen scrambling. Otherwise no scrambling is done.
  618. Default is True.
  619. permutations : array_like, optional
  620. Permutations used for scrambling.
  621. rng : `numpy.random.Generator`, optional
  622. Pseudorandom number generator state. When `rng` is None, a new
  623. `numpy.random.Generator` is created using entropy from the
  624. operating system. Types other than `numpy.random.Generator` are
  625. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  626. workers : int, optional
  627. Number of workers to use for parallel processing. If -1 is
  628. given all CPU threads are used. Default is 1.
  629. Returns
  630. -------
  631. sequence : list (n,)
  632. Sequence of Van der Corput.
  633. References
  634. ----------
  635. .. [1] A. B. Owen. "A randomized Halton algorithm in R",
  636. :arxiv:`1706.02808`, 2017.
  637. """
  638. if base < 2:
  639. raise ValueError("'base' must be at least 2")
  640. if scramble:
  641. if permutations is None:
  642. permutations = _van_der_corput_permutations(
  643. base=base, rng=rng
  644. )
  645. else:
  646. permutations = np.asarray(permutations)
  647. permutations = permutations.astype(np.int64)
  648. return _cy_van_der_corput_scrambled(n, base, start_index,
  649. permutations, workers)
  650. else:
  651. return _cy_van_der_corput(n, base, start_index, workers)
  652. class QMCEngine(ABC):
  653. """A generic Quasi-Monte Carlo sampler class meant for subclassing.
  654. QMCEngine is a base class to construct a specific Quasi-Monte Carlo
  655. sampler. It cannot be used directly as a sampler.
  656. Parameters
  657. ----------
  658. d : int
  659. Dimension of the parameter space.
  660. optimization : {None, "random-cd", "lloyd"}, optional
  661. Whether to use an optimization scheme to improve the quality after
  662. sampling. Note that this is a post-processing step that does not
  663. guarantee that all properties of the sample will be conserved.
  664. Default is None.
  665. * ``random-cd``: random permutations of coordinates to lower the
  666. centered discrepancy. The best sample based on the centered
  667. discrepancy is constantly updated. Centered discrepancy-based
  668. sampling shows better space-filling robustness toward 2D and 3D
  669. subprojections compared to using other discrepancy measures.
  670. * ``lloyd``: Perturb samples using a modified Lloyd-Max algorithm.
  671. The process converges to equally spaced samples.
  672. .. versionadded:: 1.10.0
  673. rng : `numpy.random.Generator`, optional
  674. Pseudorandom number generator state. When `rng` is None, a new
  675. `numpy.random.Generator` is created using entropy from the
  676. operating system. Types other than `numpy.random.Generator` are
  677. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  678. .. versionchanged:: 1.15.0
  679. As part of the `SPEC-007 <https://scientific-python.org/specs/spec-0007/>`_
  680. transition from use of `numpy.random.RandomState` to
  681. `numpy.random.Generator`, this keyword was changed from `seed` to
  682. `rng`. For an interim period, both keywords will continue to work, although
  683. only one may be specified at a time. After the interim period, function
  684. calls using the `seed` keyword will emit warnings. Following a
  685. deprecation period, the `seed` keyword will be removed.
  686. Notes
  687. -----
  688. By convention samples are distributed over the half-open interval
  689. ``[0, 1)``. Instances of the class can access the attributes: ``d`` for
  690. the dimension; and ``rng`` for the random number generator.
  691. **Subclassing**
  692. When subclassing `QMCEngine` to create a new sampler, ``__init__`` and
  693. ``random`` must be redefined.
  694. * ``__init__(d, rng=None)``: at least fix the dimension. If the sampler
  695. does not take advantage of a ``rng`` (deterministic methods like
  696. Halton), this parameter can be omitted.
  697. * ``_random(n, *, workers=1)``: draw ``n`` from the engine. ``workers``
  698. is used for parallelism. See `Halton` for example.
  699. Optionally, two other methods can be overwritten by subclasses:
  700. * ``reset``: Reset the engine to its original state.
  701. * ``fast_forward``: If the sequence is deterministic (like Halton
  702. sequence), then ``fast_forward(n)`` is skipping the ``n`` first draw.
  703. Examples
  704. --------
  705. To create a random sampler based on ``np.random.random``, we would do the
  706. following:
  707. >>> from scipy.stats import qmc
  708. >>> class RandomEngine(qmc.QMCEngine):
  709. ... def __init__(self, d, rng=None):
  710. ... super().__init__(d=d, rng=rng)
  711. ...
  712. ...
  713. ... def _random(self, n=1, *, workers=1):
  714. ... return self.rng.random((n, self.d))
  715. ...
  716. ...
  717. ... def reset(self):
  718. ... super().__init__(d=self.d, rng=self.rng_seed)
  719. ... return self
  720. ...
  721. ...
  722. ... def fast_forward(self, n):
  723. ... self.random(n)
  724. ... return self
  725. After subclassing `QMCEngine` to define the sampling strategy we want to
  726. use, we can create an instance to sample from.
  727. >>> engine = RandomEngine(2)
  728. >>> engine.random(5)
  729. array([[0.22733602, 0.31675834], # random
  730. [0.79736546, 0.67625467],
  731. [0.39110955, 0.33281393],
  732. [0.59830875, 0.18673419],
  733. [0.67275604, 0.94180287]])
  734. We can also reset the state of the generator and resample again.
  735. >>> _ = engine.reset()
  736. >>> engine.random(5)
  737. array([[0.22733602, 0.31675834], # random
  738. [0.79736546, 0.67625467],
  739. [0.39110955, 0.33281393],
  740. [0.59830875, 0.18673419],
  741. [0.67275604, 0.94180287]])
  742. """
  743. @abstractmethod
  744. @_transition_to_rng('seed', replace_doc=False)
  745. def __init__(
  746. self,
  747. d: IntNumber,
  748. *,
  749. optimization: Literal["random-cd", "lloyd"] | None = None,
  750. rng: SeedType = None
  751. ) -> None:
  752. self._initialize(d, optimization=optimization, rng=rng)
  753. # During SPEC 7 transition:
  754. # `__init__` has to be wrapped with @_transition_to_rng decorator
  755. # because it is public. Subclasses previously called `__init__`
  756. # directly, but this was problematic because arguments passed to
  757. # subclass `__init__` as `seed` would get passed to superclass
  758. # `__init__` as `rng`, rejecting `RandomState` arguments.
  759. def _initialize(
  760. self,
  761. d: IntNumber,
  762. *,
  763. optimization: Literal["random-cd", "lloyd"] | None = None,
  764. rng: SeedType = None
  765. ) -> None:
  766. if not np.issubdtype(type(d), np.integer) or d < 0:
  767. raise ValueError('d must be a non-negative integer value')
  768. self.d = d
  769. if isinstance(rng, np.random.Generator):
  770. # Spawn a Generator that we can own and reset.
  771. self.rng = _rng_spawn(rng, 1)[0]
  772. else:
  773. # Create our instance of Generator, does not need spawning
  774. # Also catch RandomState which cannot be spawned
  775. self.rng = check_random_state(rng)
  776. self.rng_seed = copy.deepcopy(self.rng)
  777. self.num_generated = 0
  778. config = {
  779. # random-cd
  780. "n_nochange": 100,
  781. "n_iters": 10_000,
  782. "rng": self.rng,
  783. # lloyd
  784. "tol": 1e-5,
  785. "maxiter": 10,
  786. "qhull_options": None,
  787. }
  788. self._optimization = optimization
  789. self.optimization_method = _select_optimizer(optimization, config)
  790. @abstractmethod
  791. def _random(
  792. self, n: IntNumber = 1, *, workers: IntNumber = 1
  793. ) -> np.ndarray:
  794. ...
  795. def random(
  796. self, n: IntNumber = 1, *, workers: IntNumber = 1
  797. ) -> np.ndarray:
  798. """Draw `n` in the half-open interval ``[0, 1)``.
  799. Parameters
  800. ----------
  801. n : int, optional
  802. Number of samples to generate in the parameter space.
  803. Default is 1.
  804. workers : int, optional
  805. Only supported with `Halton`.
  806. Number of workers to use for parallel processing. If -1 is
  807. given all CPU threads are used. Default is 1. It becomes faster
  808. than one worker for `n` greater than :math:`10^3`.
  809. Returns
  810. -------
  811. sample : array_like (n, d)
  812. QMC sample.
  813. """
  814. sample = self._random(n, workers=workers)
  815. if self.optimization_method is not None:
  816. sample = self.optimization_method(sample)
  817. self.num_generated += n
  818. return sample
  819. def integers(
  820. self,
  821. l_bounds: "npt.ArrayLike",
  822. *,
  823. u_bounds: "npt.ArrayLike | None" = None,
  824. n: IntNumber = 1,
  825. endpoint: bool = False,
  826. workers: IntNumber = 1
  827. ) -> np.ndarray:
  828. r"""
  829. Draw `n` integers from `l_bounds` (inclusive) to `u_bounds`
  830. (exclusive), or if endpoint=True, `l_bounds` (inclusive) to
  831. `u_bounds` (inclusive).
  832. Parameters
  833. ----------
  834. l_bounds : int or array-like of ints
  835. Lowest (signed) integers to be drawn (unless ``u_bounds=None``,
  836. in which case this parameter is 0 and this value is used for
  837. `u_bounds`).
  838. u_bounds : int or array-like of ints, optional
  839. If provided, one above the largest (signed) integer to be drawn
  840. (see above for behavior if ``u_bounds=None``).
  841. If array-like, must contain integer values.
  842. n : int, optional
  843. Number of samples to generate in the parameter space.
  844. Default is 1.
  845. endpoint : bool, optional
  846. If true, sample from the interval ``[l_bounds, u_bounds]`` instead
  847. of the default ``[l_bounds, u_bounds)``. Defaults is False.
  848. workers : int, optional
  849. Number of workers to use for parallel processing. If -1 is
  850. given all CPU threads are used. Only supported when using `Halton`
  851. Default is 1.
  852. Returns
  853. -------
  854. sample : array_like (n, d)
  855. QMC sample.
  856. Notes
  857. -----
  858. It is safe to just use the same ``[0, 1)`` to integer mapping
  859. with QMC that you would use with MC. You still get unbiasedness,
  860. a strong law of large numbers, an asymptotically infinite variance
  861. reduction and a finite sample variance bound.
  862. To convert a sample from :math:`[0, 1)` to :math:`[a, b), b>a`,
  863. with :math:`a` the lower bounds and :math:`b` the upper bounds,
  864. the following transformation is used:
  865. .. math::
  866. \text{floor}((b - a) \cdot \text{sample} + a)
  867. """
  868. if u_bounds is None:
  869. u_bounds = l_bounds
  870. l_bounds = 0
  871. u_bounds = np.atleast_1d(u_bounds)
  872. l_bounds = np.atleast_1d(l_bounds)
  873. if endpoint:
  874. u_bounds = u_bounds + 1
  875. if (not np.issubdtype(l_bounds.dtype, np.integer) or
  876. not np.issubdtype(u_bounds.dtype, np.integer)):
  877. message = ("'u_bounds' and 'l_bounds' must be integers or"
  878. " array-like of integers")
  879. raise ValueError(message)
  880. if isinstance(self, Halton):
  881. sample = self.random(n=n, workers=workers)
  882. else:
  883. sample = self.random(n=n)
  884. sample = scale(sample, l_bounds=l_bounds, u_bounds=u_bounds)
  885. sample = np.floor(sample).astype(np.int64)
  886. return sample
  887. def reset(self) -> "QMCEngine":
  888. """Reset the engine to base state.
  889. Returns
  890. -------
  891. engine : QMCEngine
  892. Engine reset to its base state.
  893. """
  894. rng = copy.deepcopy(self.rng_seed)
  895. self.rng = check_random_state(rng)
  896. self.num_generated = 0
  897. return self
  898. def fast_forward(self, n: IntNumber) -> "QMCEngine":
  899. """Fast-forward the sequence by `n` positions.
  900. Parameters
  901. ----------
  902. n : int
  903. Number of points to skip in the sequence.
  904. Returns
  905. -------
  906. engine : QMCEngine
  907. Engine reset to its base state.
  908. """
  909. self.random(n=n)
  910. return self
  911. class Halton(QMCEngine):
  912. """Halton sequence.
  913. Pseudo-random number generator that generalize the Van der Corput sequence
  914. for multiple dimensions. The Halton sequence uses the base-two Van der
  915. Corput sequence for the first dimension, base-three for its second and
  916. base-:math:`p` for its :math:`n`-dimension, with :math:`p` the
  917. :math:`n`'th prime.
  918. Parameters
  919. ----------
  920. d : int
  921. Dimension of the parameter space.
  922. scramble : bool, optional
  923. If True, use random scrambling from [2]_. Otherwise no scrambling
  924. is done.
  925. Default is True.
  926. optimization : {None, "random-cd", "lloyd"}, optional
  927. Whether to use an optimization scheme to improve the quality after
  928. sampling. Note that this is a post-processing step that does not
  929. guarantee that all properties of the sample will be conserved.
  930. Default is None.
  931. * ``random-cd``: random permutations of coordinates to lower the
  932. centered discrepancy. The best sample based on the centered
  933. discrepancy is constantly updated. Centered discrepancy-based
  934. sampling shows better space-filling robustness toward 2D and 3D
  935. subprojections compared to using other discrepancy measures.
  936. * ``lloyd``: Perturb samples using a modified Lloyd-Max algorithm.
  937. The process converges to equally spaced samples.
  938. .. versionadded:: 1.10.0
  939. rng : `numpy.random.Generator`, optional
  940. Pseudorandom number generator state. When `rng` is None, a new
  941. `numpy.random.Generator` is created using entropy from the
  942. operating system. Types other than `numpy.random.Generator` are
  943. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  944. .. versionchanged:: 1.15.0
  945. As part of the `SPEC-007 <https://scientific-python.org/specs/spec-0007/>`_
  946. transition from use of `numpy.random.RandomState` to
  947. `numpy.random.Generator`, this keyword was changed from `seed` to
  948. `rng`. For an interim period, both keywords will continue to work, although
  949. only one may be specified at a time. After the interim period, function
  950. calls using the `seed` keyword will emit warnings. Following a
  951. deprecation period, the `seed` keyword will be removed.
  952. Notes
  953. -----
  954. The Halton sequence has severe striping artifacts for even modestly
  955. large dimensions. These can be ameliorated by scrambling. Scrambling
  956. also supports replication-based error estimates and extends
  957. applicability to unbounded integrands.
  958. References
  959. ----------
  960. .. [1] Halton, "On the efficiency of certain quasi-random sequences of
  961. points in evaluating multi-dimensional integrals", Numerische
  962. Mathematik, 1960.
  963. .. [2] A. B. Owen. "A randomized Halton algorithm in R",
  964. :arxiv:`1706.02808`, 2017.
  965. Examples
  966. --------
  967. Generate samples from a low discrepancy sequence of Halton.
  968. >>> from scipy.stats import qmc
  969. >>> sampler = qmc.Halton(d=2, scramble=False)
  970. >>> sample = sampler.random(n=5)
  971. >>> sample
  972. array([[0. , 0. ],
  973. [0.5 , 0.33333333],
  974. [0.25 , 0.66666667],
  975. [0.75 , 0.11111111],
  976. [0.125 , 0.44444444]])
  977. Compute the quality of the sample using the discrepancy criterion.
  978. >>> qmc.discrepancy(sample)
  979. 0.088893711419753
  980. If some wants to continue an existing design, extra points can be obtained
  981. by calling again `random`. Alternatively, you can skip some points like:
  982. >>> _ = sampler.fast_forward(5)
  983. >>> sample_continued = sampler.random(n=5)
  984. >>> sample_continued
  985. array([[0.3125 , 0.37037037],
  986. [0.8125 , 0.7037037 ],
  987. [0.1875 , 0.14814815],
  988. [0.6875 , 0.48148148],
  989. [0.4375 , 0.81481481]])
  990. Finally, samples can be scaled to bounds.
  991. >>> l_bounds = [0, 2]
  992. >>> u_bounds = [10, 5]
  993. >>> qmc.scale(sample_continued, l_bounds, u_bounds)
  994. array([[3.125 , 3.11111111],
  995. [8.125 , 4.11111111],
  996. [1.875 , 2.44444444],
  997. [6.875 , 3.44444444],
  998. [4.375 , 4.44444444]])
  999. """
  1000. @_transition_to_rng('seed', replace_doc=False)
  1001. def __init__(
  1002. self, d: IntNumber, *, scramble: bool = True,
  1003. optimization: Literal["random-cd", "lloyd"] | None = None,
  1004. rng: SeedType = None
  1005. ) -> None:
  1006. # Used in `scipy.integrate.qmc_quad`
  1007. self._init_quad = {'d': d, 'scramble': True,
  1008. 'optimization': optimization}
  1009. super()._initialize(d=d, optimization=optimization, rng=rng)
  1010. # important to have ``type(bdim) == int`` for performance reason
  1011. self.base = [int(bdim) for bdim in n_primes(d)]
  1012. self.scramble = scramble
  1013. self._initialize_permutations()
  1014. def _initialize_permutations(self) -> None:
  1015. """Initialize permutations for all Van der Corput sequences.
  1016. Permutations are only needed for scrambling.
  1017. """
  1018. self._permutations: list = [None] * len(self.base)
  1019. if self.scramble:
  1020. for i, bdim in enumerate(self.base):
  1021. permutations = _van_der_corput_permutations(
  1022. base=bdim, rng=self.rng
  1023. )
  1024. self._permutations[i] = permutations
  1025. def _random(
  1026. self, n: IntNumber = 1, *, workers: IntNumber = 1
  1027. ) -> np.ndarray:
  1028. """Draw `n` in the half-open interval ``[0, 1)``.
  1029. Parameters
  1030. ----------
  1031. n : int, optional
  1032. Number of samples to generate in the parameter space. Default is 1.
  1033. workers : int, optional
  1034. Number of workers to use for parallel processing. If -1 is
  1035. given all CPU threads are used. Default is 1. It becomes faster
  1036. than one worker for `n` greater than :math:`10^3`.
  1037. Returns
  1038. -------
  1039. sample : array_like (n, d)
  1040. QMC sample.
  1041. """
  1042. workers = _validate_workers(workers)
  1043. # Generate a sample using a Van der Corput sequence per dimension.
  1044. sample = [van_der_corput(n, bdim, start_index=self.num_generated,
  1045. scramble=self.scramble,
  1046. permutations=self._permutations[i],
  1047. workers=workers)
  1048. for i, bdim in enumerate(self.base)]
  1049. return np.array(sample).T.reshape(n, self.d)
  1050. class LatinHypercube(QMCEngine):
  1051. r"""Latin hypercube sampling (LHS).
  1052. A Latin hypercube sample [1]_ generates :math:`n` points in
  1053. :math:`[0,1)^{d}`. Each univariate marginal distribution is stratified,
  1054. placing exactly one point in :math:`[j/n, (j+1)/n)` for
  1055. :math:`j=0,1,...,n-1`. They are still applicable when :math:`n << d`.
  1056. Parameters
  1057. ----------
  1058. d : int
  1059. Dimension of the parameter space.
  1060. scramble : bool, optional
  1061. When False, center samples within cells of a multi-dimensional grid.
  1062. Otherwise, samples are randomly placed within cells of the grid.
  1063. .. note::
  1064. Setting ``scramble=False`` does not ensure deterministic output.
  1065. For that, use the `rng` parameter.
  1066. Default is True.
  1067. .. versionadded:: 1.10.0
  1068. optimization : {None, "random-cd", "lloyd"}, optional
  1069. Whether to use an optimization scheme to improve the quality after
  1070. sampling. Note that this is a post-processing step that does not
  1071. guarantee that all properties of the sample will be conserved.
  1072. Default is None.
  1073. * ``random-cd``: random permutations of coordinates to lower the
  1074. centered discrepancy. The best sample based on the centered
  1075. discrepancy is constantly updated. Centered discrepancy-based
  1076. sampling shows better space-filling robustness toward 2D and 3D
  1077. subprojections compared to using other discrepancy measures.
  1078. * ``lloyd``: Perturb samples using a modified Lloyd-Max algorithm.
  1079. The process converges to equally spaced samples.
  1080. .. versionadded:: 1.8.0
  1081. .. versionchanged:: 1.10.0
  1082. Add ``lloyd``.
  1083. strength : {1, 2}, optional
  1084. Strength of the LHS. ``strength=1`` produces a plain LHS while
  1085. ``strength=2`` produces an orthogonal array based LHS of strength 2
  1086. [7]_, [8]_. In that case, only ``n=p**2`` points can be sampled,
  1087. with ``p`` a prime number. It also constrains ``d <= p + 1``.
  1088. Default is 1.
  1089. .. versionadded:: 1.8.0
  1090. rng : `numpy.random.Generator`, optional
  1091. Pseudorandom number generator state. When `rng` is None, a new
  1092. `numpy.random.Generator` is created using entropy from the
  1093. operating system. Types other than `numpy.random.Generator` are
  1094. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  1095. .. versionchanged:: 1.15.0
  1096. As part of the `SPEC-007 <https://scientific-python.org/specs/spec-0007/>`_
  1097. transition from use of `numpy.random.RandomState` to
  1098. `numpy.random.Generator`, this keyword was changed from `seed` to
  1099. `rng`. For an interim period, both keywords will continue to work, although
  1100. only one may be specified at a time. After the interim period, function
  1101. calls using the `seed` keyword will emit warnings. Following a
  1102. deprecation period, the `seed` keyword will be removed.
  1103. See Also
  1104. --------
  1105. :ref:`quasi-monte-carlo`
  1106. Notes
  1107. -----
  1108. When LHS is used for integrating a function :math:`f` over :math:`n`,
  1109. LHS is extremely effective on integrands that are nearly additive [2]_.
  1110. With a LHS of :math:`n` points, the variance of the integral is always
  1111. lower than plain MC on :math:`n-1` points [3]_. There is a central limit
  1112. theorem for LHS on the mean and variance of the integral [4]_, but not
  1113. necessarily for optimized LHS due to the randomization.
  1114. :math:`A` is called an orthogonal array of strength :math:`t` if in each
  1115. n-row-by-t-column submatrix of :math:`A`: all :math:`p^t` possible
  1116. distinct rows occur the same number of times. The elements of :math:`A`
  1117. are in the set :math:`\{0, 1, ..., p-1\}`, also called symbols.
  1118. The constraint that :math:`p` must be a prime number is to allow modular
  1119. arithmetic. Increasing strength adds some symmetry to the sub-projections
  1120. of a sample. With strength 2, samples are symmetric along the diagonals of
  1121. 2D sub-projections. This may be undesirable, but on the other hand, the
  1122. sample dispersion is improved.
  1123. Strength 1 (plain LHS) brings an advantage over strength 0 (MC) and
  1124. strength 2 is a useful increment over strength 1. Going to strength 3 is
  1125. a smaller increment and scrambled QMC like Sobol', Halton are more
  1126. performant [7]_.
  1127. To create a LHS of strength 2, the orthogonal array :math:`A` is
  1128. randomized by applying a random, bijective map of the set of symbols onto
  1129. itself. For example, in column 0, all 0s might become 2; in column 1,
  1130. all 0s might become 1, etc.
  1131. Then, for each column :math:`i` and symbol :math:`j`, we add a plain,
  1132. one-dimensional LHS of size :math:`p` to the subarray where
  1133. :math:`A^i = j`. The resulting matrix is finally divided by :math:`p`.
  1134. References
  1135. ----------
  1136. .. [1] Mckay et al., "A Comparison of Three Methods for Selecting Values
  1137. of Input Variables in the Analysis of Output from a Computer Code."
  1138. Technometrics, 1979.
  1139. .. [2] M. Stein, "Large sample properties of simulations using Latin
  1140. hypercube sampling." Technometrics 29, no. 2: 143-151, 1987.
  1141. .. [3] A. B. Owen, "Monte Carlo variance of scrambled net quadrature."
  1142. SIAM Journal on Numerical Analysis 34, no. 5: 1884-1910, 1997
  1143. .. [4] Loh, W.-L. "On Latin hypercube sampling." The annals of statistics
  1144. 24, no. 5: 2058-2080, 1996.
  1145. .. [5] Fang et al. "Design and modeling for computer experiments".
  1146. Computer Science and Data Analysis Series, 2006.
  1147. .. [6] Damblin et al., "Numerical studies of space filling designs:
  1148. optimization of Latin Hypercube Samples and subprojection properties."
  1149. Journal of Simulation, 2013.
  1150. .. [7] A. B. Owen , "Orthogonal arrays for computer experiments,
  1151. integration and visualization." Statistica Sinica, 1992.
  1152. .. [8] B. Tang, "Orthogonal Array-Based Latin Hypercubes."
  1153. Journal of the American Statistical Association, 1993.
  1154. .. [9] Seaholm, Susan K. et al. (1988). Latin hypercube sampling and the
  1155. sensitivity analysis of a Monte Carlo epidemic model. Int J Biomed
  1156. Comput, 23(1-2), 97-112. :doi:`10.1016/0020-7101(88)90067-0`
  1157. Examples
  1158. --------
  1159. Generate samples from a Latin hypercube generator.
  1160. >>> from scipy.stats import qmc
  1161. >>> sampler = qmc.LatinHypercube(d=2)
  1162. >>> sample = sampler.random(n=5)
  1163. >>> sample
  1164. array([[0.1545328 , 0.53664833], # random
  1165. [0.84052691, 0.06474907],
  1166. [0.52177809, 0.93343721],
  1167. [0.68033825, 0.36265316],
  1168. [0.26544879, 0.61163943]])
  1169. Compute the quality of the sample using the discrepancy criterion.
  1170. >>> qmc.discrepancy(sample)
  1171. 0.0196... # random
  1172. Samples can be scaled to bounds.
  1173. >>> l_bounds = [0, 2]
  1174. >>> u_bounds = [10, 5]
  1175. >>> qmc.scale(sample, l_bounds, u_bounds)
  1176. array([[1.54532796, 3.609945 ], # random
  1177. [8.40526909, 2.1942472 ],
  1178. [5.2177809 , 4.80031164],
  1179. [6.80338249, 3.08795949],
  1180. [2.65448791, 3.83491828]])
  1181. Below are other examples showing alternative ways to construct LHS with
  1182. even better coverage of the space.
  1183. Using a base LHS as a baseline.
  1184. >>> sampler = qmc.LatinHypercube(d=2)
  1185. >>> sample = sampler.random(n=5)
  1186. >>> qmc.discrepancy(sample)
  1187. 0.0196... # random
  1188. Use the `optimization` keyword argument to produce a LHS with
  1189. lower discrepancy at higher computational cost.
  1190. >>> sampler = qmc.LatinHypercube(d=2, optimization="random-cd")
  1191. >>> sample = sampler.random(n=5)
  1192. >>> qmc.discrepancy(sample)
  1193. 0.0176... # random
  1194. Use the `strength` keyword argument to produce an orthogonal array based
  1195. LHS of strength 2. In this case, the number of sample points must be the
  1196. square of a prime number.
  1197. >>> sampler = qmc.LatinHypercube(d=2, strength=2)
  1198. >>> sample = sampler.random(n=9)
  1199. >>> qmc.discrepancy(sample)
  1200. 0.00526... # random
  1201. Options could be combined to produce an optimized centered
  1202. orthogonal array based LHS. After optimization, the result would not
  1203. be guaranteed to be of strength 2.
  1204. **Real-world example**
  1205. In [9]_, a Latin Hypercube sampling (LHS) strategy was used to sample a
  1206. parameter space to study the importance of each parameter of an epidemic
  1207. model. Such analysis is also called a sensitivity analysis.
  1208. Since the dimensionality of the problem is high (6), it is computationally
  1209. expensive to cover the space. When numerical experiments are costly, QMC
  1210. enables analysis that may not be possible if using a grid.
  1211. The six parameters of the model represented the probability of illness,
  1212. the probability of withdrawal, and four contact probabilities. The
  1213. authors assumed uniform distributions for all parameters and generated
  1214. 50 samples.
  1215. Using `scipy.stats.qmc.LatinHypercube` to replicate the protocol,
  1216. the first step is to create a sample in the unit hypercube:
  1217. >>> from scipy.stats import qmc
  1218. >>> sampler = qmc.LatinHypercube(d=6)
  1219. >>> sample = sampler.random(n=50)
  1220. Then the sample can be scaled to the appropriate bounds:
  1221. >>> l_bounds = [0.000125, 0.01, 0.0025, 0.05, 0.47, 0.7]
  1222. >>> u_bounds = [0.000375, 0.03, 0.0075, 0.15, 0.87, 0.9]
  1223. >>> sample_scaled = qmc.scale(sample, l_bounds, u_bounds)
  1224. Such a sample was used to run the model 50 times, and a polynomial
  1225. response surface was constructed. This allowed the authors to study the
  1226. relative importance of each parameter across the range of possibilities
  1227. of every other parameter.
  1228. In this computer experiment, they showed a 14-fold reduction in the
  1229. number of samples required to maintain an error below 2% on their
  1230. response surface when compared to a grid sampling.
  1231. """
  1232. @_transition_to_rng('seed', replace_doc=False)
  1233. def __init__(
  1234. self, d: IntNumber, *,
  1235. scramble: bool = True,
  1236. strength: int = 1,
  1237. optimization: Literal["random-cd", "lloyd"] | None = None,
  1238. rng: SeedType = None
  1239. ) -> None:
  1240. # Used in `scipy.integrate.qmc_quad`
  1241. self._init_quad = {'d': d, 'scramble': True, 'strength': strength,
  1242. 'optimization': optimization}
  1243. super()._initialize(d=d, rng=rng, optimization=optimization)
  1244. self.scramble = scramble
  1245. lhs_method_strength = {
  1246. 1: self._random_lhs,
  1247. 2: self._random_oa_lhs
  1248. }
  1249. try:
  1250. self.lhs_method: Callable = lhs_method_strength[strength]
  1251. except KeyError as exc:
  1252. message = (f"{strength!r} is not a valid strength. It must be one"
  1253. f" of {set(lhs_method_strength)!r}")
  1254. raise ValueError(message) from exc
  1255. def _random(
  1256. self, n: IntNumber = 1, *, workers: IntNumber = 1
  1257. ) -> np.ndarray:
  1258. lhs = self.lhs_method(n)
  1259. return lhs
  1260. def _random_lhs(self, n: IntNumber = 1) -> np.ndarray:
  1261. """Base LHS algorithm."""
  1262. if not self.scramble:
  1263. samples: np.ndarray | float = 0.5
  1264. else:
  1265. samples = self.rng.uniform(size=(n, self.d))
  1266. perms = np.tile(np.arange(1, n + 1),
  1267. (self.d, 1)) # type: ignore[arg-type]
  1268. for i in range(self.d):
  1269. self.rng.shuffle(perms[i, :])
  1270. perms = perms.T
  1271. samples = (perms - samples) / n
  1272. return samples
  1273. def _random_oa_lhs(self, n: IntNumber = 4) -> np.ndarray:
  1274. """Orthogonal array based LHS of strength 2."""
  1275. p = np.sqrt(n).astype(int)
  1276. n_row = p**2
  1277. n_col = p + 1
  1278. primes = primes_from_2_to(p + 1)
  1279. if p not in primes or n != n_row:
  1280. raise ValueError(
  1281. "n is not the square of a prime number. Close"
  1282. f" values are {primes[-2:]**2}"
  1283. )
  1284. if self.d > p + 1:
  1285. raise ValueError("n is too small for d. Must be n > (d-1)**2")
  1286. oa_sample = np.zeros(shape=(n_row, n_col), dtype=int)
  1287. # OA of strength 2
  1288. arrays = np.tile(np.arange(p), (2, 1))
  1289. oa_sample[:, :2] = np.stack(np.meshgrid(*arrays),
  1290. axis=-1).reshape(-1, 2)
  1291. for p_ in range(1, p):
  1292. oa_sample[:, 2+p_-1] = np.mod(oa_sample[:, 0]
  1293. + p_*oa_sample[:, 1], p)
  1294. # scramble the OA
  1295. oa_sample_ = np.empty(shape=(n_row, n_col), dtype=int)
  1296. for j in range(n_col):
  1297. perms = self.rng.permutation(p)
  1298. oa_sample_[:, j] = perms[oa_sample[:, j]]
  1299. oa_sample = oa_sample_
  1300. # following is making a scrambled OA into an OA-LHS
  1301. oa_lhs_sample = np.zeros(shape=(n_row, n_col))
  1302. lhs_engine = LatinHypercube(d=1, scramble=self.scramble, strength=1,
  1303. rng=self.rng) # type: QMCEngine
  1304. for j in range(n_col):
  1305. for k in range(p):
  1306. idx = oa_sample[:, j] == k
  1307. lhs = lhs_engine.random(p).flatten()
  1308. oa_lhs_sample[:, j][idx] = lhs + oa_sample[:, j][idx]
  1309. oa_lhs_sample /= p
  1310. return oa_lhs_sample[:, :self.d]
  1311. class Sobol(QMCEngine):
  1312. """Engine for generating (scrambled) Sobol' sequences.
  1313. Sobol' sequences are low-discrepancy, quasi-random numbers. Points
  1314. can be drawn using two methods:
  1315. * `random_base2`: safely draw :math:`n=2^m` points. This method
  1316. guarantees the balance properties of the sequence.
  1317. * `random`: draw an arbitrary number of points from the
  1318. sequence. See warning below.
  1319. Parameters
  1320. ----------
  1321. d : int
  1322. Dimensionality of the sequence. Max dimensionality is 21201.
  1323. scramble : bool, optional
  1324. If True, use LMS+shift scrambling. Otherwise, no scrambling is done.
  1325. Default is True.
  1326. bits : int, optional
  1327. Number of bits of the generator. Control the maximum number of points
  1328. that can be generated, which is ``2**bits``. Maximal value is 64.
  1329. It does not correspond to the return type, which is always
  1330. ``np.float64`` to prevent points from repeating themselves.
  1331. Default is None, which for backward compatibility, corresponds to 30.
  1332. .. versionadded:: 1.9.0
  1333. optimization : {None, "random-cd", "lloyd"}, optional
  1334. Whether to use an optimization scheme to improve the quality after
  1335. sampling. Note that this is a post-processing step that does not
  1336. guarantee that all properties of the sample will be conserved.
  1337. Default is None.
  1338. * ``random-cd``: random permutations of coordinates to lower the
  1339. centered discrepancy. The best sample based on the centered
  1340. discrepancy is constantly updated. Centered discrepancy-based
  1341. sampling shows better space-filling robustness toward 2D and 3D
  1342. subprojections compared to using other discrepancy measures.
  1343. * ``lloyd``: Perturb samples using a modified Lloyd-Max algorithm.
  1344. The process converges to equally spaced samples.
  1345. .. versionadded:: 1.10.0
  1346. rng : `numpy.random.Generator`, optional
  1347. Pseudorandom number generator state. When `rng` is None, a new
  1348. `numpy.random.Generator` is created using entropy from the
  1349. operating system. Types other than `numpy.random.Generator` are
  1350. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  1351. .. versionchanged:: 1.15.0
  1352. As part of the `SPEC-007 <https://scientific-python.org/specs/spec-0007/>`_
  1353. transition from use of `numpy.random.RandomState` to
  1354. `numpy.random.Generator`, this keyword was changed from `seed` to
  1355. `rng`. For an interim period, both keywords will continue to work, although
  1356. only one may be specified at a time. After the interim period, function
  1357. calls using the `seed` keyword will emit warnings. Following a
  1358. deprecation period, the `seed` keyword will be removed.
  1359. Notes
  1360. -----
  1361. Sobol' sequences [1]_ provide :math:`n=2^m` low discrepancy points in
  1362. :math:`[0,1)^{d}`. Scrambling them [3]_ makes them suitable for singular
  1363. integrands, provides a means of error estimation, and can improve their
  1364. rate of convergence. The scrambling strategy which is implemented is a
  1365. (left) linear matrix scramble (LMS) followed by a digital random shift
  1366. (LMS+shift) [2]_.
  1367. There are many versions of Sobol' sequences depending on their
  1368. 'direction numbers'. This code uses direction numbers from [4]_. Hence,
  1369. the maximum number of dimension is 21201. The direction numbers have been
  1370. precomputed with search criterion 6 and can be retrieved at
  1371. https://web.maths.unsw.edu.au/~fkuo/sobol/.
  1372. .. warning::
  1373. Sobol' sequences are a quadrature rule and they lose their balance
  1374. properties if one uses a sample size that is not a power of 2, or skips
  1375. the first point, or thins the sequence [5]_.
  1376. If :math:`n=2^m` points are not enough then one should take :math:`2^M`
  1377. points for :math:`M>m`. When scrambling, the number R of independent
  1378. replicates does not have to be a power of 2.
  1379. Sobol' sequences are generated to some number :math:`B` of bits.
  1380. After :math:`2^B` points have been generated, the sequence would
  1381. repeat. Hence, an error is raised.
  1382. The number of bits can be controlled with the parameter `bits`.
  1383. References
  1384. ----------
  1385. .. [1] I. M. Sobol', "The distribution of points in a cube and the accurate
  1386. evaluation of integrals." Zh. Vychisl. Mat. i Mat. Phys., 7:784-802,
  1387. 1967.
  1388. .. [2] J. Matousek, "On the L2-discrepancy for anchored boxes."
  1389. J. of Complexity 14, 527-556, 1998.
  1390. .. [3] Art B. Owen, "Scrambling Sobol and Niederreiter-Xing points."
  1391. Journal of Complexity, 14(4):466-489, December 1998.
  1392. .. [4] S. Joe and F. Y. Kuo, "Constructing sobol sequences with better
  1393. two-dimensional projections." SIAM Journal on Scientific Computing,
  1394. 30(5):2635-2654, 2008.
  1395. .. [5] Art B. Owen, "On dropping the first Sobol' point."
  1396. :arxiv:`2008.08051`, 2020.
  1397. Examples
  1398. --------
  1399. Generate samples from a low discrepancy sequence of Sobol'.
  1400. >>> from scipy.stats import qmc
  1401. >>> sampler = qmc.Sobol(d=2, scramble=False)
  1402. >>> sample = sampler.random_base2(m=3)
  1403. >>> sample
  1404. array([[0. , 0. ],
  1405. [0.5 , 0.5 ],
  1406. [0.75 , 0.25 ],
  1407. [0.25 , 0.75 ],
  1408. [0.375, 0.375],
  1409. [0.875, 0.875],
  1410. [0.625, 0.125],
  1411. [0.125, 0.625]])
  1412. Compute the quality of the sample using the discrepancy criterion.
  1413. >>> qmc.discrepancy(sample)
  1414. 0.013882107204860938
  1415. To continue an existing design, extra points can be obtained
  1416. by calling again `random_base2`. Alternatively, you can skip some
  1417. points like:
  1418. >>> _ = sampler.reset()
  1419. >>> _ = sampler.fast_forward(4)
  1420. >>> sample_continued = sampler.random_base2(m=2)
  1421. >>> sample_continued
  1422. array([[0.375, 0.375],
  1423. [0.875, 0.875],
  1424. [0.625, 0.125],
  1425. [0.125, 0.625]])
  1426. Finally, samples can be scaled to bounds.
  1427. >>> l_bounds = [0, 2]
  1428. >>> u_bounds = [10, 5]
  1429. >>> qmc.scale(sample_continued, l_bounds, u_bounds)
  1430. array([[3.75 , 3.125],
  1431. [8.75 , 4.625],
  1432. [6.25 , 2.375],
  1433. [1.25 , 3.875]])
  1434. """
  1435. MAXDIM: ClassVar[int] = _MAXDIM
  1436. @_transition_to_rng('seed', replace_doc=False)
  1437. def __init__(
  1438. self, d: IntNumber, *, scramble: bool = True,
  1439. bits: IntNumber | None = None, rng: SeedType = None,
  1440. optimization: Literal["random-cd", "lloyd"] | None = None
  1441. ) -> None:
  1442. # Used in `scipy.integrate.qmc_quad`
  1443. self._init_quad = {'d': d, 'scramble': True, 'bits': bits,
  1444. 'optimization': optimization}
  1445. super()._initialize(d=d, optimization=optimization, rng=rng)
  1446. if d > self.MAXDIM:
  1447. raise ValueError(
  1448. f"Maximum supported dimensionality is {self.MAXDIM}."
  1449. )
  1450. self.bits = bits
  1451. self.dtype_i: type
  1452. self.scramble = scramble
  1453. if self.bits is None:
  1454. self.bits = 30
  1455. if self.bits <= 32:
  1456. self.dtype_i = np.uint32
  1457. elif 32 < self.bits <= 64:
  1458. self.dtype_i = np.uint64
  1459. else:
  1460. raise ValueError("Maximum supported 'bits' is 64")
  1461. self.maxn = 2**self.bits
  1462. # v is d x maxbit matrix
  1463. self._sv: np.ndarray = np.zeros((d, self.bits), dtype=self.dtype_i)
  1464. _initialize_v(self._sv, dim=d, bits=self.bits)
  1465. if not scramble:
  1466. self._shift: np.ndarray = np.zeros(d, dtype=self.dtype_i)
  1467. else:
  1468. # scramble self._shift and self._sv
  1469. self._scramble()
  1470. self._quasi = self._shift.copy()
  1471. # normalization constant with the largest possible number
  1472. # calculate in Python to not overflow int with 2**64
  1473. self._scale = 1.0 / 2 ** self.bits
  1474. self._first_point = (self._quasi * self._scale).reshape(1, -1)
  1475. # explicit casting to float64
  1476. self._first_point = self._first_point.astype(np.float64)
  1477. def _scramble(self) -> None:
  1478. """Scramble the sequence using LMS+shift."""
  1479. # Generate shift vector
  1480. self._shift = np.dot(
  1481. rng_integers(self.rng, 2, size=(self.d, self.bits),
  1482. dtype=self.dtype_i),
  1483. 2 ** np.arange(self.bits, dtype=self.dtype_i),
  1484. )
  1485. # Generate lower triangular matrices (stacked across dimensions)
  1486. ltm = np.tril(rng_integers(self.rng, 2,
  1487. size=(self.d, self.bits, self.bits),
  1488. dtype=self.dtype_i))
  1489. _cscramble(
  1490. dim=self.d, bits=self.bits, # type: ignore[arg-type]
  1491. ltm=ltm, sv=self._sv
  1492. )
  1493. def _random(
  1494. self, n: IntNumber = 1, *, workers: IntNumber = 1
  1495. ) -> np.ndarray:
  1496. """Draw next point(s) in the Sobol' sequence.
  1497. Parameters
  1498. ----------
  1499. n : int, optional
  1500. Number of samples to generate in the parameter space. Default is 1.
  1501. Returns
  1502. -------
  1503. sample : array_like (n, d)
  1504. Sobol' sample.
  1505. """
  1506. sample: np.ndarray = np.empty((n, self.d), dtype=np.float64)
  1507. if n == 0:
  1508. return sample
  1509. total_n = self.num_generated + n
  1510. if total_n > self.maxn:
  1511. msg = (
  1512. f"At most 2**{self.bits}={self.maxn} distinct points can be "
  1513. f"generated. {self.num_generated} points have been previously "
  1514. f"generated, then: n={self.num_generated}+{n}={total_n}. "
  1515. )
  1516. if self.bits != 64:
  1517. msg += "Consider increasing `bits`."
  1518. raise ValueError(msg)
  1519. if self.num_generated == 0:
  1520. # verify n is 2**n
  1521. if not (n & (n - 1) == 0):
  1522. warnings.warn("The balance properties of Sobol' points require"
  1523. " n to be a power of 2.", stacklevel=3)
  1524. if n == 1:
  1525. sample = self._first_point
  1526. else:
  1527. _draw(
  1528. n=n - 1, num_gen=self.num_generated, dim=self.d,
  1529. scale=self._scale, sv=self._sv, quasi=self._quasi,
  1530. sample=sample
  1531. )
  1532. sample = np.concatenate(
  1533. [self._first_point, sample]
  1534. )[:n]
  1535. else:
  1536. _draw(
  1537. n=n, num_gen=self.num_generated - 1, dim=self.d,
  1538. scale=self._scale, sv=self._sv, quasi=self._quasi,
  1539. sample=sample
  1540. )
  1541. return sample
  1542. def random_base2(self, m: IntNumber) -> np.ndarray:
  1543. """Draw point(s) from the Sobol' sequence.
  1544. This function draws :math:`n=2^m` points in the parameter space
  1545. ensuring the balance properties of the sequence.
  1546. Parameters
  1547. ----------
  1548. m : int
  1549. Logarithm in base 2 of the number of samples; i.e., n = 2^m.
  1550. Returns
  1551. -------
  1552. sample : array_like (n, d)
  1553. Sobol' sample.
  1554. """
  1555. n = 2 ** m
  1556. total_n = self.num_generated + n
  1557. if not (total_n & (total_n - 1) == 0):
  1558. raise ValueError('The balance properties of Sobol\' points require '
  1559. f'n to be a power of 2. {self.num_generated} points '
  1560. 'have been previously generated, then: '
  1561. f'n={self.num_generated}+2**{m}={total_n}. '
  1562. 'If you still want to do this, the function '
  1563. '\'Sobol.random()\' can be used.'
  1564. )
  1565. return self.random(n)
  1566. def reset(self) -> "Sobol":
  1567. """Reset the engine to base state.
  1568. Returns
  1569. -------
  1570. engine : Sobol
  1571. Engine reset to its base state.
  1572. """
  1573. super().reset()
  1574. self._quasi = self._shift.copy()
  1575. return self
  1576. def fast_forward(self, n: IntNumber) -> "Sobol":
  1577. """Fast-forward the sequence by `n` positions.
  1578. Parameters
  1579. ----------
  1580. n : int
  1581. Number of points to skip in the sequence.
  1582. Returns
  1583. -------
  1584. engine : Sobol
  1585. The fast-forwarded engine.
  1586. """
  1587. if self.num_generated == 0:
  1588. _fast_forward(
  1589. n=n - 1, num_gen=self.num_generated, dim=self.d,
  1590. sv=self._sv, quasi=self._quasi
  1591. )
  1592. else:
  1593. _fast_forward(
  1594. n=n, num_gen=self.num_generated - 1, dim=self.d,
  1595. sv=self._sv, quasi=self._quasi
  1596. )
  1597. self.num_generated += n
  1598. return self
  1599. class PoissonDisk(QMCEngine):
  1600. """Poisson disk sampling.
  1601. Parameters
  1602. ----------
  1603. d : int
  1604. Dimension of the parameter space.
  1605. radius : float
  1606. Minimal distance to keep between points when sampling new candidates.
  1607. hypersphere : {"volume", "surface"}, optional
  1608. Sampling strategy to generate potential candidates to be added in the
  1609. final sample. Default is "volume".
  1610. * ``volume``: original Bridson algorithm as described in [1]_.
  1611. New candidates are sampled *within* the hypersphere.
  1612. * ``surface``: only sample the surface of the hypersphere.
  1613. ncandidates : int
  1614. Number of candidates to sample per iteration. More candidates result
  1615. in a denser sampling as more candidates can be accepted per iteration.
  1616. optimization : {None, "random-cd", "lloyd"}, optional
  1617. Whether to use an optimization scheme to improve the quality after
  1618. sampling. Note that this is a post-processing step that does not
  1619. guarantee that all properties of the sample will be conserved.
  1620. Default is None.
  1621. * ``random-cd``: random permutations of coordinates to lower the
  1622. centered discrepancy. The best sample based on the centered
  1623. discrepancy is constantly updated. Centered discrepancy-based
  1624. sampling shows better space-filling robustness toward 2D and 3D
  1625. subprojections compared to using other discrepancy measures.
  1626. * ``lloyd``: Perturb samples using a modified Lloyd-Max algorithm.
  1627. The process converges to equally spaced samples.
  1628. .. versionadded:: 1.10.0
  1629. rng : `numpy.random.Generator`, optional
  1630. Pseudorandom number generator state. When `rng` is None, a new
  1631. `numpy.random.Generator` is created using entropy from the
  1632. operating system. Types other than `numpy.random.Generator` are
  1633. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  1634. .. versionchanged:: 1.15.0
  1635. As part of the `SPEC-007 <https://scientific-python.org/specs/spec-0007/>`_
  1636. transition from use of `numpy.random.RandomState` to
  1637. `numpy.random.Generator`, this keyword was changed from `seed` to
  1638. `rng`. For an interim period, both keywords will continue to work, although
  1639. only one may be specified at a time. After the interim period, function
  1640. calls using the `seed` keyword will emit warnings. Following a
  1641. deprecation period, the `seed` keyword will be removed.
  1642. l_bounds, u_bounds : array_like (d,)
  1643. Lower and upper bounds of target sample data.
  1644. Notes
  1645. -----
  1646. Poisson disk sampling is an iterative sampling strategy. Starting from
  1647. a seed sample, `ncandidates` are sampled in the hypersphere
  1648. surrounding the seed. Candidates below a certain `radius` or outside the
  1649. domain are rejected. New samples are added in a pool of sample seed. The
  1650. process stops when the pool is empty or when the number of required
  1651. samples is reached.
  1652. The maximum number of point that a sample can contain is directly linked
  1653. to the `radius`. As the dimension of the space increases, a higher radius
  1654. spreads the points further and help overcome the curse of dimensionality.
  1655. See the :ref:`quasi monte carlo tutorial <quasi-monte-carlo>` for more
  1656. details.
  1657. .. warning::
  1658. The algorithm is more suitable for low dimensions and sampling size
  1659. due to its iterative nature and memory requirements.
  1660. Selecting a small radius with a high dimension would
  1661. mean that the space could contain more samples than using lower
  1662. dimension or a bigger radius.
  1663. Some code taken from [2]_, written consent given on 31.03.2021
  1664. by the original author, Shamis, for free use in SciPy under
  1665. the 3-clause BSD.
  1666. References
  1667. ----------
  1668. .. [1] Robert Bridson, "Fast Poisson Disk Sampling in Arbitrary
  1669. Dimensions." SIGGRAPH, 2007.
  1670. .. [2] `StackOverflow <https://stackoverflow.com/questions/66047540>`__.
  1671. Examples
  1672. --------
  1673. Generate a 2D sample using a `radius` of 0.2.
  1674. >>> import numpy as np
  1675. >>> import matplotlib.pyplot as plt
  1676. >>> from matplotlib.collections import PatchCollection
  1677. >>> from scipy.stats import qmc
  1678. >>>
  1679. >>> rng = np.random.default_rng()
  1680. >>> radius = 0.2
  1681. >>> engine = qmc.PoissonDisk(d=2, radius=radius, rng=rng)
  1682. >>> sample = engine.random(20)
  1683. Visualizing the 2D sample and showing that no points are closer than
  1684. `radius`. ``radius/2`` is used to visualize non-intersecting circles.
  1685. If two samples are exactly at `radius` from each other, then their circle
  1686. of radius ``radius/2`` will touch.
  1687. >>> fig, ax = plt.subplots()
  1688. >>> _ = ax.scatter(sample[:, 0], sample[:, 1])
  1689. >>> circles = [plt.Circle((xi, yi), radius=radius/2, fill=False)
  1690. ... for xi, yi in sample]
  1691. >>> collection = PatchCollection(circles, match_original=True)
  1692. >>> ax.add_collection(collection)
  1693. >>> _ = ax.set(aspect='equal', xlabel=r'$x_1$', ylabel=r'$x_2$',
  1694. ... xlim=[0, 1], ylim=[0, 1])
  1695. >>> plt.show()
  1696. Such visualization can be seen as circle packing: how many circle can
  1697. we put in the space. It is a np-hard problem. The method `fill_space`
  1698. can be used to add samples until no more samples can be added. This is
  1699. a hard problem and parameters may need to be adjusted manually. Beware of
  1700. the dimension: as the dimensionality increases, the number of samples
  1701. required to fill the space increases exponentially
  1702. (curse-of-dimensionality).
  1703. """
  1704. @_transition_to_rng('seed', replace_doc=False)
  1705. def __init__(
  1706. self,
  1707. d: IntNumber,
  1708. *,
  1709. radius: DecimalNumber = 0.05,
  1710. hypersphere: Literal["volume", "surface"] = "volume",
  1711. ncandidates: IntNumber = 30,
  1712. optimization: Literal["random-cd", "lloyd"] | None = None,
  1713. rng: SeedType = None,
  1714. l_bounds: "npt.ArrayLike | None" = None,
  1715. u_bounds: "npt.ArrayLike | None" = None,
  1716. ) -> None:
  1717. # Used in `scipy.integrate.qmc_quad`
  1718. self._init_quad = {'d': d, 'radius': radius,
  1719. 'hypersphere': hypersphere,
  1720. 'ncandidates': ncandidates,
  1721. 'optimization': optimization}
  1722. super()._initialize(d=d, optimization=optimization, rng=rng)
  1723. hypersphere_sample = {
  1724. "volume": self._hypersphere_volume_sample,
  1725. "surface": self._hypersphere_surface_sample
  1726. }
  1727. try:
  1728. self.hypersphere_method = hypersphere_sample[hypersphere]
  1729. except KeyError as exc:
  1730. message = (
  1731. f"{hypersphere!r} is not a valid hypersphere sampling"
  1732. f" method. It must be one of {set(hypersphere_sample)!r}")
  1733. raise ValueError(message) from exc
  1734. # size of the sphere from which the samples are drawn relative to the
  1735. # size of a disk (radius)
  1736. # for the surface sampler, all new points are almost exactly 1 radius
  1737. # away from at least one existing sample +eps to avoid rejection
  1738. self.radius_factor = 2 if hypersphere == "volume" else 1.001
  1739. self.radius = radius
  1740. self.radius_squared = self.radius**2
  1741. # sample to generate per iteration in the hypersphere around center
  1742. self.ncandidates = ncandidates
  1743. if u_bounds is None:
  1744. u_bounds = np.ones(d)
  1745. if l_bounds is None:
  1746. l_bounds = np.zeros(d)
  1747. self.l_bounds, self.u_bounds = _validate_bounds(
  1748. l_bounds=l_bounds, u_bounds=u_bounds, d=int(d)
  1749. )
  1750. with np.errstate(divide='ignore'):
  1751. self.cell_size = self.radius / np.sqrt(self.d)
  1752. self.grid_size = (
  1753. np.ceil((self.u_bounds - self.l_bounds) / self.cell_size)
  1754. ).astype(int)
  1755. self._initialize_grid_pool()
  1756. def _initialize_grid_pool(self):
  1757. """Sampling pool and sample grid."""
  1758. self.sample_pool = []
  1759. # Positions of cells
  1760. # n-dim value for each grid cell
  1761. self.sample_grid = np.empty(
  1762. np.append(self.grid_size, self.d),
  1763. dtype=np.float32
  1764. )
  1765. # Initialise empty cells with NaNs
  1766. self.sample_grid.fill(np.nan)
  1767. def _random(
  1768. self, n: IntNumber = 1, *, workers: IntNumber = 1
  1769. ) -> np.ndarray:
  1770. """Draw `n` in the interval ``[l_bounds, u_bounds]``.
  1771. Note that it can return fewer samples if the space is full.
  1772. See the note section of the class.
  1773. Parameters
  1774. ----------
  1775. n : int, optional
  1776. Number of samples to generate in the parameter space. Default is 1.
  1777. Returns
  1778. -------
  1779. sample : array_like (n, d)
  1780. QMC sample.
  1781. """
  1782. if n == 0 or self.d == 0:
  1783. return np.empty((n, self.d))
  1784. def in_limits(sample: np.ndarray) -> bool:
  1785. for i in range(self.d):
  1786. if (sample[i] > self.u_bounds[i] or sample[i] < self.l_bounds[i]):
  1787. return False
  1788. return True
  1789. def in_neighborhood(candidate: np.ndarray, n: int = 2) -> bool:
  1790. """
  1791. Check if there are samples closer than ``radius_squared`` to the
  1792. `candidate` sample.
  1793. """
  1794. indices = ((candidate - self.l_bounds) / self.cell_size).astype(int)
  1795. ind_min = np.maximum(indices - n, 0)
  1796. ind_max = np.minimum(indices + n + 1, self.grid_size)
  1797. # Check if the center cell is empty
  1798. if not np.isnan(self.sample_grid[tuple(indices)][0]):
  1799. return True
  1800. a = [slice(ind_min[i], ind_max[i]) for i in range(self.d)]
  1801. # guards against: invalid value encountered in less as we are
  1802. # comparing with nan and returns False. Which is wanted.
  1803. with np.errstate(invalid='ignore'):
  1804. if np.any(
  1805. np.sum(
  1806. np.square(candidate - self.sample_grid[tuple(a)]),
  1807. axis=self.d
  1808. ) < self.radius_squared
  1809. ):
  1810. return True
  1811. return False
  1812. def add_sample(candidate: np.ndarray) -> None:
  1813. self.sample_pool.append(candidate)
  1814. indices = ((candidate - self.l_bounds) / self.cell_size).astype(int)
  1815. self.sample_grid[tuple(indices)] = candidate
  1816. curr_sample.append(candidate)
  1817. curr_sample: list[np.ndarray] = []
  1818. if len(self.sample_pool) == 0:
  1819. # the pool is being initialized with a single random sample
  1820. add_sample(self.rng.uniform(self.l_bounds, self.u_bounds))
  1821. num_drawn = 1
  1822. else:
  1823. num_drawn = 0
  1824. # exhaust sample pool to have up to n sample
  1825. while len(self.sample_pool) and num_drawn < n:
  1826. # select a sample from the available pool
  1827. idx_center = rng_integers(self.rng, len(self.sample_pool))
  1828. center = self.sample_pool[idx_center]
  1829. del self.sample_pool[idx_center]
  1830. # generate candidates around the center sample
  1831. candidates = self.hypersphere_method(
  1832. center, self.radius * self.radius_factor, self.ncandidates
  1833. )
  1834. # keep candidates that satisfy some conditions
  1835. for candidate in candidates:
  1836. if in_limits(candidate) and not in_neighborhood(candidate):
  1837. add_sample(candidate)
  1838. num_drawn += 1
  1839. if num_drawn >= n:
  1840. break
  1841. self.num_generated += num_drawn
  1842. return np.array(curr_sample)
  1843. def fill_space(self) -> np.ndarray:
  1844. """Draw ``n`` samples in the interval ``[l_bounds, u_bounds]``.
  1845. Unlike `random`, this method will try to add points until
  1846. the space is full. Depending on ``candidates`` (and to a lesser extent
  1847. other parameters), some empty areas can still be present in the sample.
  1848. .. warning::
  1849. This can be extremely slow in high dimensions or if the
  1850. ``radius`` is very small-with respect to the dimensionality.
  1851. Returns
  1852. -------
  1853. sample : array_like (n, d)
  1854. QMC sample.
  1855. """
  1856. return self.random(np.inf) # type: ignore[arg-type]
  1857. def reset(self) -> "PoissonDisk":
  1858. """Reset the engine to base state.
  1859. Returns
  1860. -------
  1861. engine : PoissonDisk
  1862. Engine reset to its base state.
  1863. """
  1864. super().reset()
  1865. self._initialize_grid_pool()
  1866. return self
  1867. def _hypersphere_volume_sample(
  1868. self, center: np.ndarray, radius: DecimalNumber,
  1869. candidates: IntNumber = 1
  1870. ) -> np.ndarray:
  1871. """Uniform sampling within hypersphere."""
  1872. # should remove samples within r/2
  1873. x = self.rng.standard_normal(size=(candidates, self.d))
  1874. ssq = np.sum(x**2, axis=1)
  1875. fr = radius * gammainc(self.d/2, ssq/2)**(1/self.d) / np.sqrt(ssq)
  1876. fr_tiled = np.tile(
  1877. fr.reshape(-1, 1), (1, self.d) # type: ignore[arg-type]
  1878. )
  1879. p = center + np.multiply(x, fr_tiled)
  1880. return p
  1881. def _hypersphere_surface_sample(
  1882. self, center: np.ndarray, radius: DecimalNumber,
  1883. candidates: IntNumber = 1
  1884. ) -> np.ndarray:
  1885. """Uniform sampling on the hypersphere's surface."""
  1886. vec = self.rng.standard_normal(size=(candidates, self.d))
  1887. vec /= np.linalg.norm(vec, axis=1)[:, None]
  1888. p = center + np.multiply(vec, radius)
  1889. return p
  1890. class MultivariateNormalQMC:
  1891. r"""QMC sampling from a multivariate Normal :math:`N(\mu, \Sigma)`.
  1892. Parameters
  1893. ----------
  1894. mean : array_like (d,)
  1895. The mean vector. Where ``d`` is the dimension.
  1896. cov : array_like (d, d), optional
  1897. The covariance matrix. If omitted, use `cov_root` instead.
  1898. If both `cov` and `cov_root` are omitted, use the identity matrix.
  1899. cov_root : array_like (d, d'), optional
  1900. A root decomposition of the covariance matrix, where ``d'`` may be less
  1901. than ``d`` if the covariance is not full rank. If omitted, use `cov`.
  1902. inv_transform : bool, optional
  1903. If True, use inverse transform instead of Box-Muller. Default is True.
  1904. engine : QMCEngine, optional
  1905. Quasi-Monte Carlo engine sampler. If None, `Sobol` is used.
  1906. rng : `numpy.random.Generator`, optional
  1907. Pseudorandom number generator state. When `rng` is None, a new
  1908. `numpy.random.Generator` is created using entropy from the
  1909. operating system. Types other than `numpy.random.Generator` are
  1910. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  1911. .. versionchanged:: 1.15.0
  1912. As part of the `SPEC-007 <https://scientific-python.org/specs/spec-0007/>`_
  1913. transition from use of `numpy.random.RandomState` to
  1914. `numpy.random.Generator`, this keyword was changed from `seed` to
  1915. `rng`. For an interim period, both keywords will continue to work, although
  1916. only one may be specified at a time. After the interim period, function
  1917. calls using the `seed` keyword will emit warnings. Following a
  1918. deprecation period, the `seed` keyword will be removed.
  1919. Examples
  1920. --------
  1921. >>> import matplotlib.pyplot as plt
  1922. >>> from scipy.stats import qmc
  1923. >>> dist = qmc.MultivariateNormalQMC(mean=[0, 5], cov=[[1, 0], [0, 1]])
  1924. >>> sample = dist.random(512)
  1925. >>> _ = plt.scatter(sample[:, 0], sample[:, 1])
  1926. >>> plt.show()
  1927. """
  1928. @_transition_to_rng('seed', replace_doc=False)
  1929. def __init__(
  1930. self,
  1931. mean: "npt.ArrayLike",
  1932. cov: "npt.ArrayLike | None" = None,
  1933. *,
  1934. cov_root: "npt.ArrayLike | None" = None,
  1935. inv_transform: bool = True,
  1936. engine: QMCEngine | None = None,
  1937. rng: SeedType = None,
  1938. ) -> None:
  1939. mean = np.asarray(np.atleast_1d(mean))
  1940. d = mean.shape[0]
  1941. if cov is not None:
  1942. # covariance matrix provided
  1943. cov = np.asarray(np.atleast_2d(cov))
  1944. # check for square/symmetric cov matrix and mean vector has the
  1945. # same d
  1946. if not mean.shape[0] == cov.shape[0]:
  1947. raise ValueError("Dimension mismatch between mean and "
  1948. "covariance.")
  1949. if not np.allclose(cov, cov.transpose()):
  1950. raise ValueError("Covariance matrix is not symmetric.")
  1951. # compute Cholesky decomp; if it fails, do the eigen decomposition
  1952. try:
  1953. cov_root = np.linalg.cholesky(cov).transpose()
  1954. except np.linalg.LinAlgError:
  1955. eigval, eigvec = np.linalg.eigh(cov)
  1956. if not np.all(eigval >= -1.0e-8):
  1957. raise ValueError("Covariance matrix not PSD.")
  1958. eigval = np.clip(eigval, 0.0, None)
  1959. cov_root = (eigvec * np.sqrt(eigval)).transpose()
  1960. elif cov_root is not None:
  1961. # root decomposition provided
  1962. cov_root = np.atleast_2d(cov_root)
  1963. if not mean.shape[0] == cov_root.shape[0]:
  1964. raise ValueError("Dimension mismatch between mean and "
  1965. "covariance.")
  1966. else:
  1967. # corresponds to identity covariance matrix
  1968. cov_root = None
  1969. self._inv_transform = inv_transform
  1970. if not inv_transform:
  1971. # to apply Box-Muller, we need an even number of dimensions
  1972. engine_dim = 2 * math.ceil(d / 2)
  1973. else:
  1974. engine_dim = d
  1975. if engine is None:
  1976. # Need this during SPEC 7 transition to prevent `RandomState`
  1977. # from being passed via `rng`.
  1978. kwarg = "seed" if isinstance(rng, np.random.RandomState) else "rng"
  1979. kwargs = {kwarg: rng}
  1980. self.engine = Sobol(
  1981. d=engine_dim, scramble=True, bits=30, **kwargs
  1982. ) # type: QMCEngine
  1983. elif isinstance(engine, QMCEngine):
  1984. if engine.d != engine_dim:
  1985. raise ValueError("Dimension of `engine` must be consistent"
  1986. " with dimensions of mean and covariance."
  1987. " If `inv_transform` is False, it must be"
  1988. " an even number.")
  1989. self.engine = engine
  1990. else:
  1991. raise ValueError("`engine` must be an instance of "
  1992. "`scipy.stats.qmc.QMCEngine` or `None`.")
  1993. self._mean = mean
  1994. self._corr_matrix = cov_root
  1995. self._d = d
  1996. def random(self, n: IntNumber = 1) -> np.ndarray:
  1997. """Draw `n` QMC samples from the multivariate Normal.
  1998. Parameters
  1999. ----------
  2000. n : int, optional
  2001. Number of samples to generate in the parameter space. Default is 1.
  2002. Returns
  2003. -------
  2004. sample : array_like (n, d)
  2005. Sample.
  2006. """
  2007. base_samples = self._standard_normal_samples(n)
  2008. return self._correlate(base_samples)
  2009. def _correlate(self, base_samples: np.ndarray) -> np.ndarray:
  2010. if self._corr_matrix is not None:
  2011. return base_samples @ self._corr_matrix + self._mean
  2012. else:
  2013. # avoid multiplying with identity here
  2014. return base_samples + self._mean
  2015. def _standard_normal_samples(self, n: IntNumber = 1) -> np.ndarray:
  2016. """Draw `n` QMC samples from the standard Normal :math:`N(0, I_d)`.
  2017. Parameters
  2018. ----------
  2019. n : int, optional
  2020. Number of samples to generate in the parameter space. Default is 1.
  2021. Returns
  2022. -------
  2023. sample : array_like (n, d)
  2024. Sample.
  2025. """
  2026. # get base samples
  2027. samples = self.engine.random(n)
  2028. if self._inv_transform:
  2029. # apply inverse transform
  2030. # (values to close to 0/1 result in inf values)
  2031. return stats.norm.ppf(0.5 + (1 - 1e-10) * (samples - 0.5)) # type: ignore[attr-defined] # noqa: E501
  2032. else:
  2033. # apply Box-Muller transform (note: indexes starting from 1)
  2034. even = np.arange(0, samples.shape[-1], 2)
  2035. Rs = np.sqrt(-2 * np.log(samples[:, even]))
  2036. thetas = 2 * math.pi * samples[:, 1 + even]
  2037. cos = np.cos(thetas)
  2038. sin = np.sin(thetas)
  2039. transf_samples = np.stack([Rs * cos, Rs * sin],
  2040. -1).reshape(n, -1)
  2041. # make sure we only return the number of dimension requested
  2042. return transf_samples[:, : self._d]
  2043. class MultinomialQMC:
  2044. r"""QMC sampling from a multinomial distribution.
  2045. Parameters
  2046. ----------
  2047. pvals : array_like (k,)
  2048. Vector of probabilities of size ``k``, where ``k`` is the number
  2049. of categories. Elements must be non-negative and sum to 1.
  2050. n_trials : int
  2051. Number of trials.
  2052. engine : QMCEngine, optional
  2053. Quasi-Monte Carlo engine sampler. If None, `Sobol` is used.
  2054. rng : `numpy.random.Generator`, optional
  2055. Pseudorandom number generator state. When `rng` is None, a new
  2056. `numpy.random.Generator` is created using entropy from the
  2057. operating system. Types other than `numpy.random.Generator` are
  2058. passed to `numpy.random.default_rng` to instantiate a ``Generator``.
  2059. .. versionchanged:: 1.15.0
  2060. As part of the `SPEC-007 <https://scientific-python.org/specs/spec-0007/>`_
  2061. transition from use of `numpy.random.RandomState` to
  2062. `numpy.random.Generator`, this keyword was changed from `seed` to
  2063. `rng`. For an interim period, both keywords will continue to work, although
  2064. only one may be specified at a time. After the interim period, function
  2065. calls using the `seed` keyword will emit warnings. Following a
  2066. deprecation period, the `seed` keyword will be removed.
  2067. Examples
  2068. --------
  2069. Let's define 3 categories and for a given sample, the sum of the trials
  2070. of each category is 8. The number of trials per category is determined
  2071. by the `pvals` associated to each category.
  2072. Then, we sample this distribution 64 times.
  2073. >>> import matplotlib.pyplot as plt
  2074. >>> from scipy.stats import qmc
  2075. >>> dist = qmc.MultinomialQMC(
  2076. ... pvals=[0.2, 0.4, 0.4], n_trials=10, engine=qmc.Halton(d=1)
  2077. ... )
  2078. >>> sample = dist.random(64)
  2079. We can plot the sample and verify that the median of number of trials
  2080. for each category is following the `pvals`. That would be
  2081. ``pvals * n_trials = [2, 4, 4]``.
  2082. >>> fig, ax = plt.subplots()
  2083. >>> ax.yaxis.get_major_locator().set_params(integer=True)
  2084. >>> _ = ax.boxplot(sample)
  2085. >>> ax.set(xlabel="Categories", ylabel="Trials")
  2086. >>> plt.show()
  2087. """
  2088. @_transition_to_rng('seed', replace_doc=False)
  2089. def __init__(
  2090. self,
  2091. pvals: "npt.ArrayLike",
  2092. n_trials: IntNumber,
  2093. *,
  2094. engine: QMCEngine | None = None,
  2095. rng: SeedType = None,
  2096. ) -> None:
  2097. self.pvals = np.atleast_1d(np.asarray(pvals))
  2098. if np.min(pvals) < 0:
  2099. raise ValueError('Elements of pvals must be non-negative.')
  2100. if not np.isclose(np.sum(pvals), 1):
  2101. raise ValueError('Elements of pvals must sum to 1.')
  2102. self.n_trials = n_trials
  2103. if engine is None:
  2104. # Need this during SPEC 7 transition to prevent `RandomState`
  2105. # from being passed via `rng`.
  2106. kwarg = "seed" if isinstance(rng, np.random.RandomState) else "rng"
  2107. kwargs = {kwarg: rng}
  2108. self.engine = Sobol(
  2109. d=1, scramble=True, bits=30, **kwargs
  2110. ) # type: QMCEngine
  2111. elif isinstance(engine, QMCEngine):
  2112. if engine.d != 1:
  2113. raise ValueError("Dimension of `engine` must be 1.")
  2114. self.engine = engine
  2115. else:
  2116. raise ValueError("`engine` must be an instance of "
  2117. "`scipy.stats.qmc.QMCEngine` or `None`.")
  2118. def random(self, n: IntNumber = 1) -> np.ndarray:
  2119. """Draw `n` QMC samples from the multinomial distribution.
  2120. Parameters
  2121. ----------
  2122. n : int, optional
  2123. Number of samples to generate in the parameter space. Default is 1.
  2124. Returns
  2125. -------
  2126. samples : array_like (n, pvals)
  2127. Sample.
  2128. """
  2129. sample = np.empty((n, len(self.pvals)))
  2130. for i in range(n):
  2131. base_draws = self.engine.random(self.n_trials).ravel()
  2132. p_cumulative = np.empty_like(self.pvals, dtype=float)
  2133. _fill_p_cumulative(np.array(self.pvals, dtype=float), p_cumulative)
  2134. sample_ = np.zeros_like(self.pvals, dtype=np.intp)
  2135. _categorize(base_draws, p_cumulative, sample_)
  2136. sample[i] = sample_
  2137. return sample
  2138. def _select_optimizer(
  2139. optimization: Literal["random-cd", "lloyd"] | None, config: dict
  2140. ) -> Callable | None:
  2141. """A factory for optimization methods."""
  2142. optimization_method: dict[str, Callable] = {
  2143. "random-cd": _random_cd,
  2144. "lloyd": _lloyd_centroidal_voronoi_tessellation
  2145. }
  2146. optimizer: partial | None
  2147. if optimization is not None:
  2148. try:
  2149. optimization = optimization.lower() # type: ignore[assignment]
  2150. optimizer_ = optimization_method[optimization]
  2151. except KeyError as exc:
  2152. message = (f"{optimization!r} is not a valid optimization"
  2153. f" method. It must be one of"
  2154. f" {set(optimization_method)!r}")
  2155. raise ValueError(message) from exc
  2156. # config
  2157. optimizer = partial(optimizer_, **config)
  2158. else:
  2159. optimizer = None
  2160. return optimizer
  2161. def _random_cd(
  2162. best_sample: np.ndarray, n_iters: int, n_nochange: int, rng: GeneratorType,
  2163. **kwargs: dict
  2164. ) -> np.ndarray:
  2165. """Optimal LHS on CD.
  2166. Create a base LHS and do random permutations of coordinates to
  2167. lower the centered discrepancy.
  2168. Because it starts with a normal LHS, it also works with the
  2169. `scramble` keyword argument.
  2170. Two stopping criterion are used to stop the algorithm: at most,
  2171. `n_iters` iterations are performed; or if there is no improvement
  2172. for `n_nochange` consecutive iterations.
  2173. """
  2174. del kwargs # only use keywords which are defined, needed by factory
  2175. n, d = best_sample.shape
  2176. if d == 0 or n == 0:
  2177. return np.empty((n, d))
  2178. if d == 1 or n == 1:
  2179. # discrepancy measures are invariant under permuting factors and runs
  2180. return best_sample
  2181. best_disc = discrepancy(best_sample)
  2182. bounds = ([0, d - 1],
  2183. [0, n - 1],
  2184. [0, n - 1])
  2185. n_nochange_ = 0
  2186. n_iters_ = 0
  2187. while n_nochange_ < n_nochange and n_iters_ < n_iters:
  2188. n_iters_ += 1
  2189. col = rng_integers(rng, *bounds[0], endpoint=True) # type: ignore[misc]
  2190. row_1 = rng_integers(rng, *bounds[1], endpoint=True) # type: ignore[misc]
  2191. row_2 = rng_integers(rng, *bounds[2], endpoint=True) # type: ignore[misc]
  2192. disc = _perturb_discrepancy(best_sample,
  2193. row_1, row_2, col,
  2194. best_disc)
  2195. if disc < best_disc:
  2196. best_sample[row_1, col], best_sample[row_2, col] = (
  2197. best_sample[row_2, col], best_sample[row_1, col])
  2198. best_disc = disc
  2199. n_nochange_ = 0
  2200. else:
  2201. n_nochange_ += 1
  2202. return best_sample
  2203. def _l1_norm(sample: np.ndarray) -> float:
  2204. return distance.pdist(sample, 'cityblock').min()
  2205. def _lloyd_iteration(
  2206. sample: np.ndarray,
  2207. decay: float,
  2208. qhull_options: str
  2209. ) -> np.ndarray:
  2210. """Lloyd-Max algorithm iteration.
  2211. Based on the implementation of Stéfan van der Walt:
  2212. https://github.com/stefanv/lloyd
  2213. which is:
  2214. Copyright (c) 2021-04-21 Stéfan van der Walt
  2215. https://github.com/stefanv/lloyd
  2216. MIT License
  2217. Parameters
  2218. ----------
  2219. sample : array_like (n, d)
  2220. The sample to iterate on.
  2221. decay : float
  2222. Relaxation decay. A positive value would move the samples toward
  2223. their centroid, and negative value would move them away.
  2224. 1 would move the samples to their centroid.
  2225. qhull_options : str
  2226. Additional options to pass to Qhull. See Qhull manual
  2227. for details. (Default: "Qbb Qc Qz Qj Qx" for ndim > 4 and
  2228. "Qbb Qc Qz Qj" otherwise.)
  2229. Returns
  2230. -------
  2231. sample : array_like (n, d)
  2232. The sample after an iteration of Lloyd's algorithm.
  2233. """
  2234. new_sample = np.empty_like(sample)
  2235. voronoi = Voronoi(sample, qhull_options=qhull_options)
  2236. for ii, idx in enumerate(voronoi.point_region):
  2237. # the region is a series of indices into self.voronoi.vertices
  2238. # remove samples at infinity, designated by index -1
  2239. region = [i for i in voronoi.regions[idx] if i != -1]
  2240. # get the vertices for this region
  2241. verts = voronoi.vertices[region]
  2242. # clipping would be wrong, we need to intersect
  2243. # verts = np.clip(verts, 0, 1)
  2244. # move samples towards centroids:
  2245. # Centroid in n-D is the mean for uniformly distributed nodes
  2246. # of a geometry.
  2247. centroid = np.mean(verts, axis=0)
  2248. new_sample[ii] = sample[ii] + (centroid - sample[ii]) * decay
  2249. # only update sample to centroid within the region
  2250. is_valid = np.all(np.logical_and(new_sample >= 0, new_sample <= 1), axis=1)
  2251. sample[is_valid] = new_sample[is_valid]
  2252. return sample
  2253. def _lloyd_centroidal_voronoi_tessellation(
  2254. sample: "npt.ArrayLike",
  2255. *,
  2256. tol: DecimalNumber = 1e-5,
  2257. maxiter: IntNumber = 10,
  2258. qhull_options: str | None = None,
  2259. **kwargs: dict
  2260. ) -> np.ndarray:
  2261. """Approximate Centroidal Voronoi Tessellation.
  2262. Perturb samples in N-dimensions using Lloyd-Max algorithm.
  2263. Parameters
  2264. ----------
  2265. sample : array_like (n, d)
  2266. The sample to iterate on. With ``n`` the number of samples and ``d``
  2267. the dimension. Samples must be in :math:`[0, 1]^d`, with ``d>=2``.
  2268. tol : float, optional
  2269. Tolerance for termination. If the min of the L1-norm over the samples
  2270. changes less than `tol`, it stops the algorithm. Default is 1e-5.
  2271. maxiter : int, optional
  2272. Maximum number of iterations. It will stop the algorithm even if
  2273. `tol` is above the threshold.
  2274. Too many iterations tend to cluster the samples as a hypersphere.
  2275. Default is 10.
  2276. qhull_options : str, optional
  2277. Additional options to pass to Qhull. See Qhull manual
  2278. for details. (Default: "Qbb Qc Qz Qj Qx" for ndim > 4 and
  2279. "Qbb Qc Qz Qj" otherwise.)
  2280. Returns
  2281. -------
  2282. sample : array_like (n, d)
  2283. The sample after being processed by Lloyd-Max algorithm.
  2284. Notes
  2285. -----
  2286. Lloyd-Max algorithm is an iterative process with the purpose of improving
  2287. the dispersion of samples. For given sample: (i) compute a Voronoi
  2288. Tessellation; (ii) find the centroid of each Voronoi cell; (iii) move the
  2289. samples toward the centroid of their respective cell. See [1]_, [2]_.
  2290. A relaxation factor is used to control how fast samples can move at each
  2291. iteration. This factor is starting at 2 and ending at 1 after `maxiter`
  2292. following an exponential decay.
  2293. The process converges to equally spaced samples. It implies that measures
  2294. like the discrepancy could suffer from too many iterations. On the other
  2295. hand, L1 and L2 distances should improve. This is especially true with
  2296. QMC methods which tend to favor the discrepancy over other criteria.
  2297. .. note::
  2298. The current implementation does not intersect the Voronoi Tessellation
  2299. with the boundaries. This implies that for a low number of samples,
  2300. empirically below 20, no Voronoi cell is touching the boundaries.
  2301. Hence, samples cannot be moved close to the boundaries.
  2302. Further improvements could consider the samples at infinity so that
  2303. all boundaries are segments of some Voronoi cells. This would fix
  2304. the computation of the centroid position.
  2305. .. warning::
  2306. The Voronoi Tessellation step is expensive and quickly becomes
  2307. intractable with dimensions as low as 10 even for a sample
  2308. of size as low as 1000.
  2309. .. versionadded:: 1.9.0
  2310. References
  2311. ----------
  2312. .. [1] Lloyd. "Least Squares Quantization in PCM".
  2313. IEEE Transactions on Information Theory, 1982.
  2314. .. [2] Max J. "Quantizing for minimum distortion".
  2315. IEEE Transactions on Information Theory, 1960.
  2316. Examples
  2317. --------
  2318. >>> import numpy as np
  2319. >>> from scipy.spatial import distance
  2320. >>> from scipy.stats._qmc import _lloyd_centroidal_voronoi_tessellation
  2321. >>> rng = np.random.default_rng()
  2322. >>> sample = rng.random((128, 2))
  2323. .. note::
  2324. The samples need to be in :math:`[0, 1]^d`. `scipy.stats.qmc.scale`
  2325. can be used to scale the samples from their
  2326. original bounds to :math:`[0, 1]^d`. And back to their original bounds.
  2327. Compute the quality of the sample using the L1 criterion.
  2328. >>> def l1_norm(sample):
  2329. ... return distance.pdist(sample, 'cityblock').min()
  2330. >>> l1_norm(sample)
  2331. 0.00161... # random
  2332. Now process the sample using Lloyd's algorithm and check the improvement
  2333. on the L1. The value should increase.
  2334. >>> sample = _lloyd_centroidal_voronoi_tessellation(sample)
  2335. >>> l1_norm(sample)
  2336. 0.0278... # random
  2337. """
  2338. del kwargs # only use keywords which are defined, needed by factory
  2339. sample = np.asarray(sample).copy()
  2340. if not sample.ndim == 2:
  2341. raise ValueError('`sample` is not a 2D array')
  2342. if not sample.shape[1] >= 2:
  2343. raise ValueError('`sample` dimension is not >= 2')
  2344. # Checking that sample is within the hypercube
  2345. if (sample.max() > 1.) or (sample.min() < 0.):
  2346. raise ValueError('`sample` is not in unit hypercube')
  2347. if qhull_options is None:
  2348. qhull_options = 'Qbb Qc Qz QJ'
  2349. if sample.shape[1] >= 5:
  2350. qhull_options += ' Qx'
  2351. # Fit an exponential to be 2 at 0 and 1 at `maxiter`.
  2352. # The decay is used for relaxation.
  2353. # analytical solution for y=exp(-maxiter/x) - 0.1
  2354. root = -maxiter / np.log(0.1)
  2355. decay = [np.exp(-x / root)+0.9 for x in range(maxiter)]
  2356. l1_old = _l1_norm(sample=sample)
  2357. for i in range(maxiter):
  2358. sample = _lloyd_iteration(
  2359. sample=sample, decay=decay[i],
  2360. qhull_options=qhull_options,
  2361. )
  2362. l1_new = _l1_norm(sample=sample)
  2363. if abs(l1_new - l1_old) < tol:
  2364. break
  2365. else:
  2366. l1_old = l1_new
  2367. return sample
  2368. def _validate_workers(workers: IntNumber = 1) -> IntNumber:
  2369. """Validate `workers` based on platform and value.
  2370. Parameters
  2371. ----------
  2372. workers : int, optional
  2373. Number of workers to use for parallel processing. If -1 is
  2374. given all CPU threads are used. Default is 1.
  2375. Returns
  2376. -------
  2377. Workers : int
  2378. Number of CPU used by the algorithm
  2379. """
  2380. workers = int(workers)
  2381. if workers == -1:
  2382. workers = os.cpu_count() # type: ignore[assignment]
  2383. if workers is None:
  2384. raise NotImplementedError(
  2385. "Cannot determine the number of cpus using os.cpu_count(), "
  2386. "cannot use -1 for the number of workers"
  2387. )
  2388. elif workers <= 0:
  2389. raise ValueError(f"Invalid number of workers: {workers}, must be -1 "
  2390. "or > 0")
  2391. return workers
  2392. def _validate_bounds(
  2393. l_bounds: "npt.ArrayLike", u_bounds: "npt.ArrayLike", d: int
  2394. ) -> "tuple[npt.NDArray[np.generic], npt.NDArray[np.generic]]":
  2395. """Bounds input validation.
  2396. Parameters
  2397. ----------
  2398. l_bounds, u_bounds : array_like (d,)
  2399. Lower and upper bounds.
  2400. d : int
  2401. Dimension to use for broadcasting.
  2402. Returns
  2403. -------
  2404. l_bounds, u_bounds : array_like (d,)
  2405. Lower and upper bounds.
  2406. """
  2407. try:
  2408. lower = np.broadcast_to(l_bounds, d)
  2409. upper = np.broadcast_to(u_bounds, d)
  2410. except ValueError as exc:
  2411. msg = ("'l_bounds' and 'u_bounds' must be broadcastable and respect"
  2412. " the sample dimension")
  2413. raise ValueError(msg) from exc
  2414. if not np.all(lower < upper):
  2415. raise ValueError("Bounds are not consistent 'l_bounds' < 'u_bounds'")
  2416. return lower, upper