distance.py 88 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797179817991800180118021803180418051806180718081809181018111812181318141815181618171818181918201821182218231824182518261827182818291830183118321833183418351836183718381839184018411842184318441845184618471848184918501851185218531854185518561857185818591860186118621863186418651866186718681869187018711872187318741875187618771878187918801881188218831884188518861887188818891890189118921893189418951896189718981899190019011902190319041905190619071908190919101911191219131914191519161917191819191920192119221923192419251926192719281929193019311932193319341935193619371938193919401941194219431944194519461947194819491950195119521953195419551956195719581959196019611962196319641965196619671968196919701971197219731974197519761977197819791980198119821983198419851986198719881989199019911992199319941995199619971998199920002001200220032004200520062007200820092010201120122013201420152016201720182019202020212022202320242025202620272028202920302031203220332034203520362037203820392040204120422043204420452046204720482049205020512052205320542055205620572058205920602061206220632064206520662067206820692070207120722073207420752076207720782079208020812082208320842085208620872088208920902091209220932094209520962097209820992100210121022103210421052106210721082109211021112112211321142115211621172118211921202121212221232124212521262127212821292130213121322133213421352136213721382139214021412142214321442145214621472148214921502151215221532154215521562157215821592160216121622163216421652166216721682169217021712172217321742175217621772178217921802181218221832184218521862187218821892190219121922193219421952196219721982199220022012202220322042205220622072208220922102211221222132214221522162217221822192220222122222223222422252226222722282229223022312232223322342235223622372238223922402241224222432244224522462247224822492250225122522253225422552256225722582259226022612262226322642265226622672268226922702271227222732274227522762277227822792280228122822283228422852286228722882289229022912292229322942295229622972298229923002301230223032304230523062307230823092310231123122313231423152316231723182319232023212322232323242325232623272328232923302331233223332334233523362337233823392340234123422343234423452346234723482349235023512352235323542355235623572358235923602361236223632364236523662367236823692370237123722373237423752376237723782379238023812382238323842385238623872388238923902391239223932394239523962397239823992400240124022403240424052406240724082409241024112412241324142415241624172418241924202421242224232424242524262427242824292430243124322433243424352436243724382439244024412442244324442445244624472448244924502451245224532454245524562457245824592460246124622463246424652466246724682469247024712472247324742475247624772478247924802481248224832484248524862487248824892490249124922493249424952496249724982499250025012502250325042505250625072508250925102511251225132514251525162517251825192520252125222523252425252526252725282529253025312532253325342535253625372538253925402541254225432544254525462547254825492550255125522553255425552556255725582559256025612562256325642565256625672568256925702571257225732574257525762577257825792580258125822583258425852586258725882589259025912592259325942595259625972598259926002601260226032604260526062607260826092610261126122613261426152616261726182619262026212622262326242625262626272628262926302631263226332634263526362637263826392640264126422643264426452646264726482649265026512652265326542655265626572658265926602661266226632664266526662667266826692670267126722673267426752676267726782679268026812682268326842685268626872688268926902691269226932694269526962697269826992700270127022703270427052706270727082709271027112712271327142715271627172718271927202721272227232724272527262727272827292730273127322733273427352736273727382739274027412742274327442745274627472748274927502751275227532754275527562757275827592760276127622763276427652766276727682769277027712772277327742775277627772778277927802781278227832784278527862787278827892790279127922793279427952796279727982799280028012802280328042805280628072808280928102811281228132814281528162817281828192820282128222823282428252826282728282829283028312832283328342835283628372838283928402841284228432844284528462847284828492850285128522853285428552856285728582859286028612862286328642865286628672868286928702871287228732874287528762877287828792880288128822883288428852886288728882889289028912892289328942895289628972898289929002901290229032904290529062907290829092910291129122913291429152916291729182919292029212922292329242925292629272928292929302931293229332934293529362937293829392940294129422943294429452946294729482949295029512952
  1. """
  2. Distance computations (:mod:`scipy.spatial.distance`)
  3. =====================================================
  4. .. sectionauthor:: Damian Eads
  5. Function reference
  6. ------------------
  7. Distance matrix computation from a collection of raw observation vectors
  8. stored in a rectangular array.
  9. .. autosummary::
  10. :toctree: generated/
  11. pdist -- pairwise distances between observation vectors.
  12. cdist -- distances between two collections of observation vectors
  13. squareform -- convert distance matrix to a condensed one and vice versa
  14. directed_hausdorff -- directed Hausdorff distance between arrays
  15. Predicates for checking the validity of distance matrices, both
  16. condensed and redundant. Also contained in this module are functions
  17. for computing the number of observations in a distance matrix.
  18. .. autosummary::
  19. :toctree: generated/
  20. is_valid_dm -- checks for a valid distance matrix
  21. is_valid_y -- checks for a valid condensed distance matrix
  22. num_obs_dm -- # of observations in a distance matrix
  23. num_obs_y -- # of observations in a condensed distance matrix
  24. Distance functions between two numeric vectors ``u`` and ``v``. Computing
  25. distances over a large collection of vectors is inefficient for these
  26. functions. Use ``pdist`` for this purpose.
  27. .. autosummary::
  28. :toctree: generated/
  29. braycurtis -- the Bray-Curtis distance.
  30. canberra -- the Canberra distance.
  31. chebyshev -- the Chebyshev distance.
  32. cityblock -- the Manhattan distance.
  33. correlation -- the Correlation distance.
  34. cosine -- the Cosine distance.
  35. euclidean -- the Euclidean distance.
  36. jensenshannon -- the Jensen-Shannon distance.
  37. mahalanobis -- the Mahalanobis distance.
  38. minkowski -- the Minkowski distance.
  39. seuclidean -- the normalized Euclidean distance.
  40. sqeuclidean -- the squared Euclidean distance.
  41. Distance functions between two boolean vectors (representing sets) ``u`` and
  42. ``v``. As in the case of numerical vectors, ``pdist`` is more efficient for
  43. computing the distances between all pairs.
  44. .. autosummary::
  45. :toctree: generated/
  46. dice -- the Dice dissimilarity.
  47. hamming -- the Hamming distance.
  48. jaccard -- the Jaccard distance.
  49. kulsinski -- the Kulsinski distance.
  50. kulczynski1 -- the Kulczynski 1 distance.
  51. rogerstanimoto -- the Rogers-Tanimoto dissimilarity.
  52. russellrao -- the Russell-Rao dissimilarity.
  53. sokalmichener -- the Sokal-Michener dissimilarity.
  54. sokalsneath -- the Sokal-Sneath dissimilarity.
  55. yule -- the Yule dissimilarity.
  56. :func:`hamming` also operates over discrete numerical vectors.
  57. """
  58. # Copyright (C) Damian Eads, 2007-2008. New BSD License.
  59. __all__ = [
  60. 'braycurtis',
  61. 'canberra',
  62. 'cdist',
  63. 'chebyshev',
  64. 'cityblock',
  65. 'correlation',
  66. 'cosine',
  67. 'dice',
  68. 'directed_hausdorff',
  69. 'euclidean',
  70. 'hamming',
  71. 'is_valid_dm',
  72. 'is_valid_y',
  73. 'jaccard',
  74. 'jensenshannon',
  75. 'kulsinski',
  76. 'kulczynski1',
  77. 'mahalanobis',
  78. 'minkowski',
  79. 'num_obs_dm',
  80. 'num_obs_y',
  81. 'pdist',
  82. 'rogerstanimoto',
  83. 'russellrao',
  84. 'seuclidean',
  85. 'sokalmichener',
  86. 'sokalsneath',
  87. 'sqeuclidean',
  88. 'squareform',
  89. 'yule'
  90. ]
  91. import warnings
  92. import numpy as np
  93. import dataclasses
  94. from typing import List, Optional, Set, Callable
  95. from functools import partial
  96. from scipy._lib._util import _asarray_validated
  97. from . import _distance_wrap
  98. from . import _hausdorff
  99. from ..linalg import norm
  100. from ..special import rel_entr
  101. from . import _distance_pybind
  102. from .._lib.deprecation import _deprecated
  103. def _copy_array_if_base_present(a):
  104. """Copy the array if its base points to a parent array."""
  105. if a.base is not None:
  106. return a.copy()
  107. return a
  108. def _correlation_cdist_wrap(XA, XB, dm, **kwargs):
  109. XA = XA - XA.mean(axis=1, keepdims=True)
  110. XB = XB - XB.mean(axis=1, keepdims=True)
  111. _distance_wrap.cdist_cosine_double_wrap(XA, XB, dm, **kwargs)
  112. def _correlation_pdist_wrap(X, dm, **kwargs):
  113. X2 = X - X.mean(axis=1, keepdims=True)
  114. _distance_wrap.pdist_cosine_double_wrap(X2, dm, **kwargs)
  115. def _convert_to_type(X, out_type):
  116. return np.ascontiguousarray(X, dtype=out_type)
  117. def _nbool_correspond_all(u, v, w=None):
  118. if u.dtype == v.dtype == bool and w is None:
  119. not_u = ~u
  120. not_v = ~v
  121. nff = (not_u & not_v).sum()
  122. nft = (not_u & v).sum()
  123. ntf = (u & not_v).sum()
  124. ntt = (u & v).sum()
  125. else:
  126. dtype = np.result_type(int, u.dtype, v.dtype)
  127. u = u.astype(dtype)
  128. v = v.astype(dtype)
  129. not_u = 1.0 - u
  130. not_v = 1.0 - v
  131. if w is not None:
  132. not_u = w * not_u
  133. u = w * u
  134. nff = (not_u * not_v).sum()
  135. nft = (not_u * v).sum()
  136. ntf = (u * not_v).sum()
  137. ntt = (u * v).sum()
  138. return (nff, nft, ntf, ntt)
  139. def _nbool_correspond_ft_tf(u, v, w=None):
  140. if u.dtype == v.dtype == bool and w is None:
  141. not_u = ~u
  142. not_v = ~v
  143. nft = (not_u & v).sum()
  144. ntf = (u & not_v).sum()
  145. else:
  146. dtype = np.result_type(int, u.dtype, v.dtype)
  147. u = u.astype(dtype)
  148. v = v.astype(dtype)
  149. not_u = 1.0 - u
  150. not_v = 1.0 - v
  151. if w is not None:
  152. not_u = w * not_u
  153. u = w * u
  154. nft = (not_u * v).sum()
  155. ntf = (u * not_v).sum()
  156. return (nft, ntf)
  157. def _validate_cdist_input(XA, XB, mA, mB, n, metric_info, **kwargs):
  158. # get supported types
  159. types = metric_info.types
  160. # choose best type
  161. typ = types[types.index(XA.dtype)] if XA.dtype in types else types[0]
  162. # validate data
  163. XA = _convert_to_type(XA, out_type=typ)
  164. XB = _convert_to_type(XB, out_type=typ)
  165. # validate kwargs
  166. _validate_kwargs = metric_info.validator
  167. if _validate_kwargs:
  168. kwargs = _validate_kwargs((XA, XB), mA + mB, n, **kwargs)
  169. return XA, XB, typ, kwargs
  170. def _validate_weight_with_size(X, m, n, **kwargs):
  171. w = kwargs.pop('w', None)
  172. if w is None:
  173. return kwargs
  174. if w.ndim != 1 or w.shape[0] != n:
  175. raise ValueError("Weights must have same size as input vector. "
  176. f"{w.shape[0]} vs. {n}")
  177. kwargs['w'] = _validate_weights(w)
  178. return kwargs
  179. def _validate_hamming_kwargs(X, m, n, **kwargs):
  180. w = kwargs.get('w', np.ones((n,), dtype='double'))
  181. if w.ndim != 1 or w.shape[0] != n:
  182. raise ValueError("Weights must have same size as input vector. %d vs. %d" % (w.shape[0], n))
  183. kwargs['w'] = _validate_weights(w)
  184. return kwargs
  185. def _validate_mahalanobis_kwargs(X, m, n, **kwargs):
  186. VI = kwargs.pop('VI', None)
  187. if VI is None:
  188. if m <= n:
  189. # There are fewer observations than the dimension of
  190. # the observations.
  191. raise ValueError("The number of observations (%d) is too "
  192. "small; the covariance matrix is "
  193. "singular. For observations with %d "
  194. "dimensions, at least %d observations "
  195. "are required." % (m, n, n + 1))
  196. if isinstance(X, tuple):
  197. X = np.vstack(X)
  198. CV = np.atleast_2d(np.cov(X.astype(np.double, copy=False).T))
  199. VI = np.linalg.inv(CV).T.copy()
  200. kwargs["VI"] = _convert_to_double(VI)
  201. return kwargs
  202. def _validate_minkowski_kwargs(X, m, n, **kwargs):
  203. kwargs = _validate_weight_with_size(X, m, n, **kwargs)
  204. if 'p' not in kwargs:
  205. kwargs['p'] = 2.
  206. else:
  207. if kwargs['p'] <= 0:
  208. raise ValueError("p must be greater than 0")
  209. return kwargs
  210. def _validate_pdist_input(X, m, n, metric_info, **kwargs):
  211. # get supported types
  212. types = metric_info.types
  213. # choose best type
  214. typ = types[types.index(X.dtype)] if X.dtype in types else types[0]
  215. # validate data
  216. X = _convert_to_type(X, out_type=typ)
  217. # validate kwargs
  218. _validate_kwargs = metric_info.validator
  219. if _validate_kwargs:
  220. kwargs = _validate_kwargs(X, m, n, **kwargs)
  221. return X, typ, kwargs
  222. def _validate_seuclidean_kwargs(X, m, n, **kwargs):
  223. V = kwargs.pop('V', None)
  224. if V is None:
  225. if isinstance(X, tuple):
  226. X = np.vstack(X)
  227. V = np.var(X.astype(np.double, copy=False), axis=0, ddof=1)
  228. else:
  229. V = np.asarray(V, order='c')
  230. if len(V.shape) != 1:
  231. raise ValueError('Variance vector V must '
  232. 'be one-dimensional.')
  233. if V.shape[0] != n:
  234. raise ValueError('Variance vector V must be of the same '
  235. 'dimension as the vectors on which the distances '
  236. 'are computed.')
  237. kwargs['V'] = _convert_to_double(V)
  238. return kwargs
  239. def _validate_vector(u, dtype=None):
  240. # XXX Is order='c' really necessary?
  241. u = np.asarray(u, dtype=dtype, order='c')
  242. if u.ndim == 1:
  243. return u
  244. raise ValueError("Input vector should be 1-D.")
  245. def _validate_weights(w, dtype=np.double):
  246. w = _validate_vector(w, dtype=dtype)
  247. if np.any(w < 0):
  248. raise ValueError("Input weights should be all non-negative")
  249. return w
  250. def directed_hausdorff(u, v, seed=0):
  251. """
  252. Compute the directed Hausdorff distance between two 2-D arrays.
  253. Distances between pairs are calculated using a Euclidean metric.
  254. Parameters
  255. ----------
  256. u : (M,N) array_like
  257. Input array.
  258. v : (O,N) array_like
  259. Input array.
  260. seed : int or None
  261. Local `numpy.random.RandomState` seed. Default is 0, a random
  262. shuffling of u and v that guarantees reproducibility.
  263. Returns
  264. -------
  265. d : double
  266. The directed Hausdorff distance between arrays `u` and `v`,
  267. index_1 : int
  268. index of point contributing to Hausdorff pair in `u`
  269. index_2 : int
  270. index of point contributing to Hausdorff pair in `v`
  271. Raises
  272. ------
  273. ValueError
  274. An exception is thrown if `u` and `v` do not have
  275. the same number of columns.
  276. Notes
  277. -----
  278. Uses the early break technique and the random sampling approach
  279. described by [1]_. Although worst-case performance is ``O(m * o)``
  280. (as with the brute force algorithm), this is unlikely in practice
  281. as the input data would have to require the algorithm to explore
  282. every single point interaction, and after the algorithm shuffles
  283. the input points at that. The best case performance is O(m), which
  284. is satisfied by selecting an inner loop distance that is less than
  285. cmax and leads to an early break as often as possible. The authors
  286. have formally shown that the average runtime is closer to O(m).
  287. .. versionadded:: 0.19.0
  288. References
  289. ----------
  290. .. [1] A. A. Taha and A. Hanbury, "An efficient algorithm for
  291. calculating the exact Hausdorff distance." IEEE Transactions On
  292. Pattern Analysis And Machine Intelligence, vol. 37 pp. 2153-63,
  293. 2015.
  294. See Also
  295. --------
  296. scipy.spatial.procrustes : Another similarity test for two data sets
  297. Examples
  298. --------
  299. Find the directed Hausdorff distance between two 2-D arrays of
  300. coordinates:
  301. >>> from scipy.spatial.distance import directed_hausdorff
  302. >>> import numpy as np
  303. >>> u = np.array([(1.0, 0.0),
  304. ... (0.0, 1.0),
  305. ... (-1.0, 0.0),
  306. ... (0.0, -1.0)])
  307. >>> v = np.array([(2.0, 0.0),
  308. ... (0.0, 2.0),
  309. ... (-2.0, 0.0),
  310. ... (0.0, -4.0)])
  311. >>> directed_hausdorff(u, v)[0]
  312. 2.23606797749979
  313. >>> directed_hausdorff(v, u)[0]
  314. 3.0
  315. Find the general (symmetric) Hausdorff distance between two 2-D
  316. arrays of coordinates:
  317. >>> max(directed_hausdorff(u, v)[0], directed_hausdorff(v, u)[0])
  318. 3.0
  319. Find the indices of the points that generate the Hausdorff distance
  320. (the Hausdorff pair):
  321. >>> directed_hausdorff(v, u)[1:]
  322. (3, 3)
  323. """
  324. u = np.asarray(u, dtype=np.float64, order='c')
  325. v = np.asarray(v, dtype=np.float64, order='c')
  326. if u.shape[1] != v.shape[1]:
  327. raise ValueError('u and v need to have the same '
  328. 'number of columns')
  329. result = _hausdorff.directed_hausdorff(u, v, seed)
  330. return result
  331. def minkowski(u, v, p=2, w=None):
  332. """
  333. Compute the Minkowski distance between two 1-D arrays.
  334. The Minkowski distance between 1-D arrays `u` and `v`,
  335. is defined as
  336. .. math::
  337. {\\|u-v\\|}_p = (\\sum{|u_i - v_i|^p})^{1/p}.
  338. \\left(\\sum{w_i(|(u_i - v_i)|^p)}\\right)^{1/p}.
  339. Parameters
  340. ----------
  341. u : (N,) array_like
  342. Input array.
  343. v : (N,) array_like
  344. Input array.
  345. p : scalar
  346. The order of the norm of the difference :math:`{\\|u-v\\|}_p`. Note
  347. that for :math:`0 < p < 1`, the triangle inequality only holds with
  348. an additional multiplicative factor, i.e. it is only a quasi-metric.
  349. w : (N,) array_like, optional
  350. The weights for each value in `u` and `v`. Default is None,
  351. which gives each value a weight of 1.0
  352. Returns
  353. -------
  354. minkowski : double
  355. The Minkowski distance between vectors `u` and `v`.
  356. Examples
  357. --------
  358. >>> from scipy.spatial import distance
  359. >>> distance.minkowski([1, 0, 0], [0, 1, 0], 1)
  360. 2.0
  361. >>> distance.minkowski([1, 0, 0], [0, 1, 0], 2)
  362. 1.4142135623730951
  363. >>> distance.minkowski([1, 0, 0], [0, 1, 0], 3)
  364. 1.2599210498948732
  365. >>> distance.minkowski([1, 1, 0], [0, 1, 0], 1)
  366. 1.0
  367. >>> distance.minkowski([1, 1, 0], [0, 1, 0], 2)
  368. 1.0
  369. >>> distance.minkowski([1, 1, 0], [0, 1, 0], 3)
  370. 1.0
  371. """
  372. u = _validate_vector(u)
  373. v = _validate_vector(v)
  374. if p <= 0:
  375. raise ValueError("p must be greater than 0")
  376. u_v = u - v
  377. if w is not None:
  378. w = _validate_weights(w)
  379. if p == 1:
  380. root_w = w
  381. elif p == 2:
  382. # better precision and speed
  383. root_w = np.sqrt(w)
  384. elif p == np.inf:
  385. root_w = (w != 0)
  386. else:
  387. root_w = np.power(w, 1/p)
  388. u_v = root_w * u_v
  389. dist = norm(u_v, ord=p)
  390. return dist
  391. def euclidean(u, v, w=None):
  392. """
  393. Computes the Euclidean distance between two 1-D arrays.
  394. The Euclidean distance between 1-D arrays `u` and `v`, is defined as
  395. .. math::
  396. {\\|u-v\\|}_2
  397. \\left(\\sum{(w_i |(u_i - v_i)|^2)}\\right)^{1/2}
  398. Parameters
  399. ----------
  400. u : (N,) array_like
  401. Input array.
  402. v : (N,) array_like
  403. Input array.
  404. w : (N,) array_like, optional
  405. The weights for each value in `u` and `v`. Default is None,
  406. which gives each value a weight of 1.0
  407. Returns
  408. -------
  409. euclidean : double
  410. The Euclidean distance between vectors `u` and `v`.
  411. Examples
  412. --------
  413. >>> from scipy.spatial import distance
  414. >>> distance.euclidean([1, 0, 0], [0, 1, 0])
  415. 1.4142135623730951
  416. >>> distance.euclidean([1, 1, 0], [0, 1, 0])
  417. 1.0
  418. """
  419. return minkowski(u, v, p=2, w=w)
  420. def sqeuclidean(u, v, w=None):
  421. """
  422. Compute the squared Euclidean distance between two 1-D arrays.
  423. The squared Euclidean distance between `u` and `v` is defined as
  424. .. math::
  425. {\\|u-v\\|}_2^2
  426. \\left(\\sum{(w_i |(u_i - v_i)|^2)}\\right)
  427. Parameters
  428. ----------
  429. u : (N,) array_like
  430. Input array.
  431. v : (N,) array_like
  432. Input array.
  433. w : (N,) array_like, optional
  434. The weights for each value in `u` and `v`. Default is None,
  435. which gives each value a weight of 1.0
  436. Returns
  437. -------
  438. sqeuclidean : double
  439. The squared Euclidean distance between vectors `u` and `v`.
  440. Examples
  441. --------
  442. >>> from scipy.spatial import distance
  443. >>> distance.sqeuclidean([1, 0, 0], [0, 1, 0])
  444. 2.0
  445. >>> distance.sqeuclidean([1, 1, 0], [0, 1, 0])
  446. 1.0
  447. """
  448. # Preserve float dtypes, but convert everything else to np.float64
  449. # for stability.
  450. utype, vtype = None, None
  451. if not (hasattr(u, "dtype") and np.issubdtype(u.dtype, np.inexact)):
  452. utype = np.float64
  453. if not (hasattr(v, "dtype") and np.issubdtype(v.dtype, np.inexact)):
  454. vtype = np.float64
  455. u = _validate_vector(u, dtype=utype)
  456. v = _validate_vector(v, dtype=vtype)
  457. u_v = u - v
  458. u_v_w = u_v # only want weights applied once
  459. if w is not None:
  460. w = _validate_weights(w)
  461. u_v_w = w * u_v
  462. return np.dot(u_v, u_v_w)
  463. def correlation(u, v, w=None, centered=True):
  464. """
  465. Compute the correlation distance between two 1-D arrays.
  466. The correlation distance between `u` and `v`, is
  467. defined as
  468. .. math::
  469. 1 - \\frac{(u - \\bar{u}) \\cdot (v - \\bar{v})}
  470. {{\\|(u - \\bar{u})\\|}_2 {\\|(v - \\bar{v})\\|}_2}
  471. where :math:`\\bar{u}` is the mean of the elements of `u`
  472. and :math:`x \\cdot y` is the dot product of :math:`x` and :math:`y`.
  473. Parameters
  474. ----------
  475. u : (N,) array_like
  476. Input array.
  477. v : (N,) array_like
  478. Input array.
  479. w : (N,) array_like, optional
  480. The weights for each value in `u` and `v`. Default is None,
  481. which gives each value a weight of 1.0
  482. centered : bool, optional
  483. If True, `u` and `v` will be centered. Default is True.
  484. Returns
  485. -------
  486. correlation : double
  487. The correlation distance between 1-D array `u` and `v`.
  488. """
  489. u = _validate_vector(u)
  490. v = _validate_vector(v)
  491. if w is not None:
  492. w = _validate_weights(w)
  493. if centered:
  494. umu = np.average(u, weights=w)
  495. vmu = np.average(v, weights=w)
  496. u = u - umu
  497. v = v - vmu
  498. uv = np.average(u * v, weights=w)
  499. uu = np.average(np.square(u), weights=w)
  500. vv = np.average(np.square(v), weights=w)
  501. dist = 1.0 - uv / np.sqrt(uu * vv)
  502. # Return absolute value to avoid small negative value due to rounding
  503. return np.abs(dist)
  504. def cosine(u, v, w=None):
  505. """
  506. Compute the Cosine distance between 1-D arrays.
  507. The Cosine distance between `u` and `v`, is defined as
  508. .. math::
  509. 1 - \\frac{u \\cdot v}
  510. {\\|u\\|_2 \\|v\\|_2}.
  511. where :math:`u \\cdot v` is the dot product of :math:`u` and
  512. :math:`v`.
  513. Parameters
  514. ----------
  515. u : (N,) array_like
  516. Input array.
  517. v : (N,) array_like
  518. Input array.
  519. w : (N,) array_like, optional
  520. The weights for each value in `u` and `v`. Default is None,
  521. which gives each value a weight of 1.0
  522. Returns
  523. -------
  524. cosine : double
  525. The Cosine distance between vectors `u` and `v`.
  526. Examples
  527. --------
  528. >>> from scipy.spatial import distance
  529. >>> distance.cosine([1, 0, 0], [0, 1, 0])
  530. 1.0
  531. >>> distance.cosine([100, 0, 0], [0, 1, 0])
  532. 1.0
  533. >>> distance.cosine([1, 1, 0], [0, 1, 0])
  534. 0.29289321881345254
  535. """
  536. # cosine distance is also referred to as 'uncentered correlation',
  537. # or 'reflective correlation'
  538. # clamp the result to 0-2
  539. return max(0, min(correlation(u, v, w=w, centered=False), 2.0))
  540. def hamming(u, v, w=None):
  541. """
  542. Compute the Hamming distance between two 1-D arrays.
  543. The Hamming distance between 1-D arrays `u` and `v`, is simply the
  544. proportion of disagreeing components in `u` and `v`. If `u` and `v` are
  545. boolean vectors, the Hamming distance is
  546. .. math::
  547. \\frac{c_{01} + c_{10}}{n}
  548. where :math:`c_{ij}` is the number of occurrences of
  549. :math:`\\mathtt{u[k]} = i` and :math:`\\mathtt{v[k]} = j` for
  550. :math:`k < n`.
  551. Parameters
  552. ----------
  553. u : (N,) array_like
  554. Input array.
  555. v : (N,) array_like
  556. Input array.
  557. w : (N,) array_like, optional
  558. The weights for each value in `u` and `v`. Default is None,
  559. which gives each value a weight of 1.0
  560. Returns
  561. -------
  562. hamming : double
  563. The Hamming distance between vectors `u` and `v`.
  564. Examples
  565. --------
  566. >>> from scipy.spatial import distance
  567. >>> distance.hamming([1, 0, 0], [0, 1, 0])
  568. 0.66666666666666663
  569. >>> distance.hamming([1, 0, 0], [1, 1, 0])
  570. 0.33333333333333331
  571. >>> distance.hamming([1, 0, 0], [2, 0, 0])
  572. 0.33333333333333331
  573. >>> distance.hamming([1, 0, 0], [3, 0, 0])
  574. 0.33333333333333331
  575. """
  576. u = _validate_vector(u)
  577. v = _validate_vector(v)
  578. if u.shape != v.shape:
  579. raise ValueError('The 1d arrays must have equal lengths.')
  580. u_ne_v = u != v
  581. if w is not None:
  582. w = _validate_weights(w)
  583. return np.average(u_ne_v, weights=w)
  584. def jaccard(u, v, w=None):
  585. """
  586. Compute the Jaccard-Needham dissimilarity between two boolean 1-D arrays.
  587. The Jaccard-Needham dissimilarity between 1-D boolean arrays `u` and `v`,
  588. is defined as
  589. .. math::
  590. \\frac{c_{TF} + c_{FT}}
  591. {c_{TT} + c_{FT} + c_{TF}}
  592. where :math:`c_{ij}` is the number of occurrences of
  593. :math:`\\mathtt{u[k]} = i` and :math:`\\mathtt{v[k]} = j` for
  594. :math:`k < n`.
  595. Parameters
  596. ----------
  597. u : (N,) array_like, bool
  598. Input array.
  599. v : (N,) array_like, bool
  600. Input array.
  601. w : (N,) array_like, optional
  602. The weights for each value in `u` and `v`. Default is None,
  603. which gives each value a weight of 1.0
  604. Returns
  605. -------
  606. jaccard : double
  607. The Jaccard distance between vectors `u` and `v`.
  608. Notes
  609. -----
  610. When both `u` and `v` lead to a `0/0` division i.e. there is no overlap
  611. between the items in the vectors the returned distance is 0. See the
  612. Wikipedia page on the Jaccard index [1]_, and this paper [2]_.
  613. .. versionchanged:: 1.2.0
  614. Previously, when `u` and `v` lead to a `0/0` division, the function
  615. would return NaN. This was changed to return 0 instead.
  616. References
  617. ----------
  618. .. [1] https://en.wikipedia.org/wiki/Jaccard_index
  619. .. [2] S. Kosub, "A note on the triangle inequality for the Jaccard
  620. distance", 2016, :arxiv:`1612.02696`
  621. Examples
  622. --------
  623. >>> from scipy.spatial import distance
  624. >>> distance.jaccard([1, 0, 0], [0, 1, 0])
  625. 1.0
  626. >>> distance.jaccard([1, 0, 0], [1, 1, 0])
  627. 0.5
  628. >>> distance.jaccard([1, 0, 0], [1, 2, 0])
  629. 0.5
  630. >>> distance.jaccard([1, 0, 0], [1, 1, 1])
  631. 0.66666666666666663
  632. """
  633. u = _validate_vector(u)
  634. v = _validate_vector(v)
  635. nonzero = np.bitwise_or(u != 0, v != 0)
  636. unequal_nonzero = np.bitwise_and((u != v), nonzero)
  637. if w is not None:
  638. w = _validate_weights(w)
  639. nonzero = w * nonzero
  640. unequal_nonzero = w * unequal_nonzero
  641. a = np.double(unequal_nonzero.sum())
  642. b = np.double(nonzero.sum())
  643. return (a / b) if b != 0 else 0
  644. @_deprecated("Kulsinski has been deprecated from scipy.spatial.distance"
  645. " in SciPy 1.9.0 and it will be removed in SciPy 1.11.0."
  646. " It is superseded by scipy.spatial.distance.kulczynski1.")
  647. def kulsinski(u, v, w=None):
  648. """
  649. Compute the Kulsinski dissimilarity between two boolean 1-D arrays.
  650. The Kulsinski dissimilarity between two boolean 1-D arrays `u` and `v`,
  651. is defined as
  652. .. math::
  653. \\frac{c_{TF} + c_{FT} - c_{TT} + n}
  654. {c_{FT} + c_{TF} + n}
  655. where :math:`c_{ij}` is the number of occurrences of
  656. :math:`\\mathtt{u[k]} = i` and :math:`\\mathtt{v[k]} = j` for
  657. :math:`k < n`.
  658. .. deprecated:: 0.12.0
  659. `kulsinski` has been deprecated from `scipy.spatial.distance` in
  660. SciPy 1.9.0 and it will be removed in SciPy 1.11.0. It is superseded
  661. by `scipy.spatial.distance.kulczynski1`.
  662. Parameters
  663. ----------
  664. u : (N,) array_like, bool
  665. Input array.
  666. v : (N,) array_like, bool
  667. Input array.
  668. w : (N,) array_like, optional
  669. The weights for each value in `u` and `v`. Default is None,
  670. which gives each value a weight of 1.0
  671. Returns
  672. -------
  673. kulsinski : double
  674. The Kulsinski distance between vectors `u` and `v`.
  675. Examples
  676. --------
  677. >>> from scipy.spatial import distance
  678. >>> distance.kulsinski([1, 0, 0], [0, 1, 0])
  679. 1.0
  680. >>> distance.kulsinski([1, 0, 0], [1, 1, 0])
  681. 0.75
  682. >>> distance.kulsinski([1, 0, 0], [2, 1, 0])
  683. 0.33333333333333331
  684. >>> distance.kulsinski([1, 0, 0], [3, 1, 0])
  685. -0.5
  686. """
  687. u = _validate_vector(u)
  688. v = _validate_vector(v)
  689. if w is None:
  690. n = float(len(u))
  691. else:
  692. w = _validate_weights(w)
  693. n = w.sum()
  694. (nff, nft, ntf, ntt) = _nbool_correspond_all(u, v, w=w)
  695. return (ntf + nft - ntt + n) / (ntf + nft + n)
  696. def kulczynski1(u, v, *, w=None):
  697. """
  698. Compute the Kulczynski 1 dissimilarity between two boolean 1-D arrays.
  699. The Kulczynski 1 dissimilarity between two boolean 1-D arrays `u` and `v`
  700. of length ``n``, is defined as
  701. .. math::
  702. \\frac{c_{11}}
  703. {c_{01} + c_{10}}
  704. where :math:`c_{ij}` is the number of occurrences of
  705. :math:`\\mathtt{u[k]} = i` and :math:`\\mathtt{v[k]} = j` for
  706. :math:`k \\in {0, 1, ..., n-1}`.
  707. Parameters
  708. ----------
  709. u : (N,) array_like, bool
  710. Input array.
  711. v : (N,) array_like, bool
  712. Input array.
  713. w : (N,) array_like, optional
  714. The weights for each value in `u` and `v`. Default is None,
  715. which gives each value a weight of 1.0
  716. Returns
  717. -------
  718. kulczynski1 : float
  719. The Kulczynski 1 distance between vectors `u` and `v`.
  720. Notes
  721. -----
  722. This measure has a minimum value of 0 and no upper limit.
  723. It is un-defined when there are no non-matches.
  724. .. versionadded:: 1.8.0
  725. References
  726. ----------
  727. .. [1] Kulczynski S. et al. Bulletin
  728. International de l'Academie Polonaise des Sciences
  729. et des Lettres, Classe des Sciences Mathematiques
  730. et Naturelles, Serie B (Sciences Naturelles). 1927;
  731. Supplement II: 57-203.
  732. Examples
  733. --------
  734. >>> from scipy.spatial import distance
  735. >>> distance.kulczynski1([1, 0, 0], [0, 1, 0])
  736. 0.0
  737. >>> distance.kulczynski1([True, False, False], [True, True, False])
  738. 1.0
  739. >>> distance.kulczynski1([True, False, False], [True])
  740. 0.5
  741. >>> distance.kulczynski1([1, 0, 0], [3, 1, 0])
  742. -3.0
  743. """
  744. u = _validate_vector(u)
  745. v = _validate_vector(v)
  746. if w is not None:
  747. w = _validate_weights(w)
  748. (_, nft, ntf, ntt) = _nbool_correspond_all(u, v, w=w)
  749. return ntt / (ntf + nft)
  750. def seuclidean(u, v, V):
  751. """
  752. Return the standardized Euclidean distance between two 1-D arrays.
  753. The standardized Euclidean distance between `u` and `v`.
  754. Parameters
  755. ----------
  756. u : (N,) array_like
  757. Input array.
  758. v : (N,) array_like
  759. Input array.
  760. V : (N,) array_like
  761. `V` is an 1-D array of component variances. It is usually computed
  762. among a larger collection vectors.
  763. Returns
  764. -------
  765. seuclidean : double
  766. The standardized Euclidean distance between vectors `u` and `v`.
  767. Examples
  768. --------
  769. >>> from scipy.spatial import distance
  770. >>> distance.seuclidean([1, 0, 0], [0, 1, 0], [0.1, 0.1, 0.1])
  771. 4.4721359549995796
  772. >>> distance.seuclidean([1, 0, 0], [0, 1, 0], [1, 0.1, 0.1])
  773. 3.3166247903553998
  774. >>> distance.seuclidean([1, 0, 0], [0, 1, 0], [10, 0.1, 0.1])
  775. 3.1780497164141406
  776. """
  777. u = _validate_vector(u)
  778. v = _validate_vector(v)
  779. V = _validate_vector(V, dtype=np.float64)
  780. if V.shape[0] != u.shape[0] or u.shape[0] != v.shape[0]:
  781. raise TypeError('V must be a 1-D array of the same dimension '
  782. 'as u and v.')
  783. return euclidean(u, v, w=1/V)
  784. def cityblock(u, v, w=None):
  785. """
  786. Compute the City Block (Manhattan) distance.
  787. Computes the Manhattan distance between two 1-D arrays `u` and `v`,
  788. which is defined as
  789. .. math::
  790. \\sum_i {\\left| u_i - v_i \\right|}.
  791. Parameters
  792. ----------
  793. u : (N,) array_like
  794. Input array.
  795. v : (N,) array_like
  796. Input array.
  797. w : (N,) array_like, optional
  798. The weights for each value in `u` and `v`. Default is None,
  799. which gives each value a weight of 1.0
  800. Returns
  801. -------
  802. cityblock : double
  803. The City Block (Manhattan) distance between vectors `u` and `v`.
  804. Examples
  805. --------
  806. >>> from scipy.spatial import distance
  807. >>> distance.cityblock([1, 0, 0], [0, 1, 0])
  808. 2
  809. >>> distance.cityblock([1, 0, 0], [0, 2, 0])
  810. 3
  811. >>> distance.cityblock([1, 0, 0], [1, 1, 0])
  812. 1
  813. """
  814. u = _validate_vector(u)
  815. v = _validate_vector(v)
  816. l1_diff = abs(u - v)
  817. if w is not None:
  818. w = _validate_weights(w)
  819. l1_diff = w * l1_diff
  820. return l1_diff.sum()
  821. def mahalanobis(u, v, VI):
  822. """
  823. Compute the Mahalanobis distance between two 1-D arrays.
  824. The Mahalanobis distance between 1-D arrays `u` and `v`, is defined as
  825. .. math::
  826. \\sqrt{ (u-v) V^{-1} (u-v)^T }
  827. where ``V`` is the covariance matrix. Note that the argument `VI`
  828. is the inverse of ``V``.
  829. Parameters
  830. ----------
  831. u : (N,) array_like
  832. Input array.
  833. v : (N,) array_like
  834. Input array.
  835. VI : array_like
  836. The inverse of the covariance matrix.
  837. Returns
  838. -------
  839. mahalanobis : double
  840. The Mahalanobis distance between vectors `u` and `v`.
  841. Examples
  842. --------
  843. >>> from scipy.spatial import distance
  844. >>> iv = [[1, 0.5, 0.5], [0.5, 1, 0.5], [0.5, 0.5, 1]]
  845. >>> distance.mahalanobis([1, 0, 0], [0, 1, 0], iv)
  846. 1.0
  847. >>> distance.mahalanobis([0, 2, 0], [0, 1, 0], iv)
  848. 1.0
  849. >>> distance.mahalanobis([2, 0, 0], [0, 1, 0], iv)
  850. 1.7320508075688772
  851. """
  852. u = _validate_vector(u)
  853. v = _validate_vector(v)
  854. VI = np.atleast_2d(VI)
  855. delta = u - v
  856. m = np.dot(np.dot(delta, VI), delta)
  857. return np.sqrt(m)
  858. def chebyshev(u, v, w=None):
  859. """
  860. Compute the Chebyshev distance.
  861. Computes the Chebyshev distance between two 1-D arrays `u` and `v`,
  862. which is defined as
  863. .. math::
  864. \\max_i {|u_i-v_i|}.
  865. Parameters
  866. ----------
  867. u : (N,) array_like
  868. Input vector.
  869. v : (N,) array_like
  870. Input vector.
  871. w : (N,) array_like, optional
  872. Unused, as 'max' is a weightless operation. Here for API consistency.
  873. Returns
  874. -------
  875. chebyshev : double
  876. The Chebyshev distance between vectors `u` and `v`.
  877. Examples
  878. --------
  879. >>> from scipy.spatial import distance
  880. >>> distance.chebyshev([1, 0, 0], [0, 1, 0])
  881. 1
  882. >>> distance.chebyshev([1, 1, 0], [0, 1, 0])
  883. 1
  884. """
  885. u = _validate_vector(u)
  886. v = _validate_vector(v)
  887. if w is not None:
  888. w = _validate_weights(w)
  889. has_weight = w > 0
  890. if has_weight.sum() < w.size:
  891. u = u[has_weight]
  892. v = v[has_weight]
  893. return max(abs(u - v))
  894. def braycurtis(u, v, w=None):
  895. """
  896. Compute the Bray-Curtis distance between two 1-D arrays.
  897. Bray-Curtis distance is defined as
  898. .. math::
  899. \\sum{|u_i-v_i|} / \\sum{|u_i+v_i|}
  900. The Bray-Curtis distance is in the range [0, 1] if all coordinates are
  901. positive, and is undefined if the inputs are of length zero.
  902. Parameters
  903. ----------
  904. u : (N,) array_like
  905. Input array.
  906. v : (N,) array_like
  907. Input array.
  908. w : (N,) array_like, optional
  909. The weights for each value in `u` and `v`. Default is None,
  910. which gives each value a weight of 1.0
  911. Returns
  912. -------
  913. braycurtis : double
  914. The Bray-Curtis distance between 1-D arrays `u` and `v`.
  915. Examples
  916. --------
  917. >>> from scipy.spatial import distance
  918. >>> distance.braycurtis([1, 0, 0], [0, 1, 0])
  919. 1.0
  920. >>> distance.braycurtis([1, 1, 0], [0, 1, 0])
  921. 0.33333333333333331
  922. """
  923. u = _validate_vector(u)
  924. v = _validate_vector(v, dtype=np.float64)
  925. l1_diff = abs(u - v)
  926. l1_sum = abs(u + v)
  927. if w is not None:
  928. w = _validate_weights(w)
  929. l1_diff = w * l1_diff
  930. l1_sum = w * l1_sum
  931. return l1_diff.sum() / l1_sum.sum()
  932. def canberra(u, v, w=None):
  933. """
  934. Compute the Canberra distance between two 1-D arrays.
  935. The Canberra distance is defined as
  936. .. math::
  937. d(u,v) = \\sum_i \\frac{|u_i-v_i|}
  938. {|u_i|+|v_i|}.
  939. Parameters
  940. ----------
  941. u : (N,) array_like
  942. Input array.
  943. v : (N,) array_like
  944. Input array.
  945. w : (N,) array_like, optional
  946. The weights for each value in `u` and `v`. Default is None,
  947. which gives each value a weight of 1.0
  948. Returns
  949. -------
  950. canberra : double
  951. The Canberra distance between vectors `u` and `v`.
  952. Notes
  953. -----
  954. When `u[i]` and `v[i]` are 0 for given i, then the fraction 0/0 = 0 is
  955. used in the calculation.
  956. Examples
  957. --------
  958. >>> from scipy.spatial import distance
  959. >>> distance.canberra([1, 0, 0], [0, 1, 0])
  960. 2.0
  961. >>> distance.canberra([1, 1, 0], [0, 1, 0])
  962. 1.0
  963. """
  964. u = _validate_vector(u)
  965. v = _validate_vector(v, dtype=np.float64)
  966. if w is not None:
  967. w = _validate_weights(w)
  968. with np.errstate(invalid='ignore'):
  969. abs_uv = abs(u - v)
  970. abs_u = abs(u)
  971. abs_v = abs(v)
  972. d = abs_uv / (abs_u + abs_v)
  973. if w is not None:
  974. d = w * d
  975. d = np.nansum(d)
  976. return d
  977. def jensenshannon(p, q, base=None, *, axis=0, keepdims=False):
  978. """
  979. Compute the Jensen-Shannon distance (metric) between
  980. two probability arrays. This is the square root
  981. of the Jensen-Shannon divergence.
  982. The Jensen-Shannon distance between two probability
  983. vectors `p` and `q` is defined as,
  984. .. math::
  985. \\sqrt{\\frac{D(p \\parallel m) + D(q \\parallel m)}{2}}
  986. where :math:`m` is the pointwise mean of :math:`p` and :math:`q`
  987. and :math:`D` is the Kullback-Leibler divergence.
  988. This routine will normalize `p` and `q` if they don't sum to 1.0.
  989. Parameters
  990. ----------
  991. p : (N,) array_like
  992. left probability vector
  993. q : (N,) array_like
  994. right probability vector
  995. base : double, optional
  996. the base of the logarithm used to compute the output
  997. if not given, then the routine uses the default base of
  998. scipy.stats.entropy.
  999. axis : int, optional
  1000. Axis along which the Jensen-Shannon distances are computed. The default
  1001. is 0.
  1002. .. versionadded:: 1.7.0
  1003. keepdims : bool, optional
  1004. If this is set to `True`, the reduced axes are left in the
  1005. result as dimensions with size one. With this option,
  1006. the result will broadcast correctly against the input array.
  1007. Default is False.
  1008. .. versionadded:: 1.7.0
  1009. Returns
  1010. -------
  1011. js : double or ndarray
  1012. The Jensen-Shannon distances between `p` and `q` along the `axis`.
  1013. Notes
  1014. -----
  1015. .. versionadded:: 1.2.0
  1016. Examples
  1017. --------
  1018. >>> from scipy.spatial import distance
  1019. >>> import numpy as np
  1020. >>> distance.jensenshannon([1.0, 0.0, 0.0], [0.0, 1.0, 0.0], 2.0)
  1021. 1.0
  1022. >>> distance.jensenshannon([1.0, 0.0], [0.5, 0.5])
  1023. 0.46450140402245893
  1024. >>> distance.jensenshannon([1.0, 0.0, 0.0], [1.0, 0.0, 0.0])
  1025. 0.0
  1026. >>> a = np.array([[1, 2, 3, 4],
  1027. ... [5, 6, 7, 8],
  1028. ... [9, 10, 11, 12]])
  1029. >>> b = np.array([[13, 14, 15, 16],
  1030. ... [17, 18, 19, 20],
  1031. ... [21, 22, 23, 24]])
  1032. >>> distance.jensenshannon(a, b, axis=0)
  1033. array([0.1954288, 0.1447697, 0.1138377, 0.0927636])
  1034. >>> distance.jensenshannon(a, b, axis=1)
  1035. array([0.1402339, 0.0399106, 0.0201815])
  1036. """
  1037. p = np.asarray(p)
  1038. q = np.asarray(q)
  1039. p = p / np.sum(p, axis=axis, keepdims=True)
  1040. q = q / np.sum(q, axis=axis, keepdims=True)
  1041. m = (p + q) / 2.0
  1042. left = rel_entr(p, m)
  1043. right = rel_entr(q, m)
  1044. left_sum = np.sum(left, axis=axis, keepdims=keepdims)
  1045. right_sum = np.sum(right, axis=axis, keepdims=keepdims)
  1046. js = left_sum + right_sum
  1047. if base is not None:
  1048. js /= np.log(base)
  1049. return np.sqrt(js / 2.0)
  1050. def yule(u, v, w=None):
  1051. """
  1052. Compute the Yule dissimilarity between two boolean 1-D arrays.
  1053. The Yule dissimilarity is defined as
  1054. .. math::
  1055. \\frac{R}{c_{TT} * c_{FF} + \\frac{R}{2}}
  1056. where :math:`c_{ij}` is the number of occurrences of
  1057. :math:`\\mathtt{u[k]} = i` and :math:`\\mathtt{v[k]} = j` for
  1058. :math:`k < n` and :math:`R = 2.0 * c_{TF} * c_{FT}`.
  1059. Parameters
  1060. ----------
  1061. u : (N,) array_like, bool
  1062. Input array.
  1063. v : (N,) array_like, bool
  1064. Input array.
  1065. w : (N,) array_like, optional
  1066. The weights for each value in `u` and `v`. Default is None,
  1067. which gives each value a weight of 1.0
  1068. Returns
  1069. -------
  1070. yule : double
  1071. The Yule dissimilarity between vectors `u` and `v`.
  1072. Examples
  1073. --------
  1074. >>> from scipy.spatial import distance
  1075. >>> distance.yule([1, 0, 0], [0, 1, 0])
  1076. 2.0
  1077. >>> distance.yule([1, 1, 0], [0, 1, 0])
  1078. 0.0
  1079. """
  1080. u = _validate_vector(u)
  1081. v = _validate_vector(v)
  1082. if w is not None:
  1083. w = _validate_weights(w)
  1084. (nff, nft, ntf, ntt) = _nbool_correspond_all(u, v, w=w)
  1085. half_R = ntf * nft
  1086. if half_R == 0:
  1087. return 0.0
  1088. else:
  1089. return float(2.0 * half_R / (ntt * nff + half_R))
  1090. def dice(u, v, w=None):
  1091. """
  1092. Compute the Dice dissimilarity between two boolean 1-D arrays.
  1093. The Dice dissimilarity between `u` and `v`, is
  1094. .. math::
  1095. \\frac{c_{TF} + c_{FT}}
  1096. {2c_{TT} + c_{FT} + c_{TF}}
  1097. where :math:`c_{ij}` is the number of occurrences of
  1098. :math:`\\mathtt{u[k]} = i` and :math:`\\mathtt{v[k]} = j` for
  1099. :math:`k < n`.
  1100. Parameters
  1101. ----------
  1102. u : (N,) array_like, bool
  1103. Input 1-D array.
  1104. v : (N,) array_like, bool
  1105. Input 1-D array.
  1106. w : (N,) array_like, optional
  1107. The weights for each value in `u` and `v`. Default is None,
  1108. which gives each value a weight of 1.0
  1109. Returns
  1110. -------
  1111. dice : double
  1112. The Dice dissimilarity between 1-D arrays `u` and `v`.
  1113. Notes
  1114. -----
  1115. This function computes the Dice dissimilarity index. To compute the
  1116. Dice similarity index, convert one to the other with similarity =
  1117. 1 - dissimilarity.
  1118. Examples
  1119. --------
  1120. >>> from scipy.spatial import distance
  1121. >>> distance.dice([1, 0, 0], [0, 1, 0])
  1122. 1.0
  1123. >>> distance.dice([1, 0, 0], [1, 1, 0])
  1124. 0.3333333333333333
  1125. >>> distance.dice([1, 0, 0], [2, 0, 0])
  1126. -0.3333333333333333
  1127. """
  1128. u = _validate_vector(u)
  1129. v = _validate_vector(v)
  1130. if w is not None:
  1131. w = _validate_weights(w)
  1132. if u.dtype == v.dtype == bool and w is None:
  1133. ntt = (u & v).sum()
  1134. else:
  1135. dtype = np.result_type(int, u.dtype, v.dtype)
  1136. u = u.astype(dtype)
  1137. v = v.astype(dtype)
  1138. if w is None:
  1139. ntt = (u * v).sum()
  1140. else:
  1141. ntt = (u * v * w).sum()
  1142. (nft, ntf) = _nbool_correspond_ft_tf(u, v, w=w)
  1143. return float((ntf + nft) / np.array(2.0 * ntt + ntf + nft))
  1144. def rogerstanimoto(u, v, w=None):
  1145. """
  1146. Compute the Rogers-Tanimoto dissimilarity between two boolean 1-D arrays.
  1147. The Rogers-Tanimoto dissimilarity between two boolean 1-D arrays
  1148. `u` and `v`, is defined as
  1149. .. math::
  1150. \\frac{R}
  1151. {c_{TT} + c_{FF} + R}
  1152. where :math:`c_{ij}` is the number of occurrences of
  1153. :math:`\\mathtt{u[k]} = i` and :math:`\\mathtt{v[k]} = j` for
  1154. :math:`k < n` and :math:`R = 2(c_{TF} + c_{FT})`.
  1155. Parameters
  1156. ----------
  1157. u : (N,) array_like, bool
  1158. Input array.
  1159. v : (N,) array_like, bool
  1160. Input array.
  1161. w : (N,) array_like, optional
  1162. The weights for each value in `u` and `v`. Default is None,
  1163. which gives each value a weight of 1.0
  1164. Returns
  1165. -------
  1166. rogerstanimoto : double
  1167. The Rogers-Tanimoto dissimilarity between vectors
  1168. `u` and `v`.
  1169. Examples
  1170. --------
  1171. >>> from scipy.spatial import distance
  1172. >>> distance.rogerstanimoto([1, 0, 0], [0, 1, 0])
  1173. 0.8
  1174. >>> distance.rogerstanimoto([1, 0, 0], [1, 1, 0])
  1175. 0.5
  1176. >>> distance.rogerstanimoto([1, 0, 0], [2, 0, 0])
  1177. -1.0
  1178. """
  1179. u = _validate_vector(u)
  1180. v = _validate_vector(v)
  1181. if w is not None:
  1182. w = _validate_weights(w)
  1183. (nff, nft, ntf, ntt) = _nbool_correspond_all(u, v, w=w)
  1184. return float(2.0 * (ntf + nft)) / float(ntt + nff + (2.0 * (ntf + nft)))
  1185. def russellrao(u, v, w=None):
  1186. """
  1187. Compute the Russell-Rao dissimilarity between two boolean 1-D arrays.
  1188. The Russell-Rao dissimilarity between two boolean 1-D arrays, `u` and
  1189. `v`, is defined as
  1190. .. math::
  1191. \\frac{n - c_{TT}}
  1192. {n}
  1193. where :math:`c_{ij}` is the number of occurrences of
  1194. :math:`\\mathtt{u[k]} = i` and :math:`\\mathtt{v[k]} = j` for
  1195. :math:`k < n`.
  1196. Parameters
  1197. ----------
  1198. u : (N,) array_like, bool
  1199. Input array.
  1200. v : (N,) array_like, bool
  1201. Input array.
  1202. w : (N,) array_like, optional
  1203. The weights for each value in `u` and `v`. Default is None,
  1204. which gives each value a weight of 1.0
  1205. Returns
  1206. -------
  1207. russellrao : double
  1208. The Russell-Rao dissimilarity between vectors `u` and `v`.
  1209. Examples
  1210. --------
  1211. >>> from scipy.spatial import distance
  1212. >>> distance.russellrao([1, 0, 0], [0, 1, 0])
  1213. 1.0
  1214. >>> distance.russellrao([1, 0, 0], [1, 1, 0])
  1215. 0.6666666666666666
  1216. >>> distance.russellrao([1, 0, 0], [2, 0, 0])
  1217. 0.3333333333333333
  1218. """
  1219. u = _validate_vector(u)
  1220. v = _validate_vector(v)
  1221. if u.dtype == v.dtype == bool and w is None:
  1222. ntt = (u & v).sum()
  1223. n = float(len(u))
  1224. elif w is None:
  1225. ntt = (u * v).sum()
  1226. n = float(len(u))
  1227. else:
  1228. w = _validate_weights(w)
  1229. ntt = (u * v * w).sum()
  1230. n = w.sum()
  1231. return float(n - ntt) / n
  1232. def sokalmichener(u, v, w=None):
  1233. """
  1234. Compute the Sokal-Michener dissimilarity between two boolean 1-D arrays.
  1235. The Sokal-Michener dissimilarity between boolean 1-D arrays `u` and `v`,
  1236. is defined as
  1237. .. math::
  1238. \\frac{R}
  1239. {S + R}
  1240. where :math:`c_{ij}` is the number of occurrences of
  1241. :math:`\\mathtt{u[k]} = i` and :math:`\\mathtt{v[k]} = j` for
  1242. :math:`k < n`, :math:`R = 2 * (c_{TF} + c_{FT})` and
  1243. :math:`S = c_{FF} + c_{TT}`.
  1244. Parameters
  1245. ----------
  1246. u : (N,) array_like, bool
  1247. Input array.
  1248. v : (N,) array_like, bool
  1249. Input array.
  1250. w : (N,) array_like, optional
  1251. The weights for each value in `u` and `v`. Default is None,
  1252. which gives each value a weight of 1.0
  1253. Returns
  1254. -------
  1255. sokalmichener : double
  1256. The Sokal-Michener dissimilarity between vectors `u` and `v`.
  1257. Examples
  1258. --------
  1259. >>> from scipy.spatial import distance
  1260. >>> distance.sokalmichener([1, 0, 0], [0, 1, 0])
  1261. 0.8
  1262. >>> distance.sokalmichener([1, 0, 0], [1, 1, 0])
  1263. 0.5
  1264. >>> distance.sokalmichener([1, 0, 0], [2, 0, 0])
  1265. -1.0
  1266. """
  1267. u = _validate_vector(u)
  1268. v = _validate_vector(v)
  1269. if w is not None:
  1270. w = _validate_weights(w)
  1271. nff, nft, ntf, ntt = _nbool_correspond_all(u, v, w=w)
  1272. return float(2.0 * (ntf + nft)) / float(ntt + nff + 2.0 * (ntf + nft))
  1273. def sokalsneath(u, v, w=None):
  1274. """
  1275. Compute the Sokal-Sneath dissimilarity between two boolean 1-D arrays.
  1276. The Sokal-Sneath dissimilarity between `u` and `v`,
  1277. .. math::
  1278. \\frac{R}
  1279. {c_{TT} + R}
  1280. where :math:`c_{ij}` is the number of occurrences of
  1281. :math:`\\mathtt{u[k]} = i` and :math:`\\mathtt{v[k]} = j` for
  1282. :math:`k < n` and :math:`R = 2(c_{TF} + c_{FT})`.
  1283. Parameters
  1284. ----------
  1285. u : (N,) array_like, bool
  1286. Input array.
  1287. v : (N,) array_like, bool
  1288. Input array.
  1289. w : (N,) array_like, optional
  1290. The weights for each value in `u` and `v`. Default is None,
  1291. which gives each value a weight of 1.0
  1292. Returns
  1293. -------
  1294. sokalsneath : double
  1295. The Sokal-Sneath dissimilarity between vectors `u` and `v`.
  1296. Examples
  1297. --------
  1298. >>> from scipy.spatial import distance
  1299. >>> distance.sokalsneath([1, 0, 0], [0, 1, 0])
  1300. 1.0
  1301. >>> distance.sokalsneath([1, 0, 0], [1, 1, 0])
  1302. 0.66666666666666663
  1303. >>> distance.sokalsneath([1, 0, 0], [2, 1, 0])
  1304. 0.0
  1305. >>> distance.sokalsneath([1, 0, 0], [3, 1, 0])
  1306. -2.0
  1307. """
  1308. u = _validate_vector(u)
  1309. v = _validate_vector(v)
  1310. if u.dtype == v.dtype == bool and w is None:
  1311. ntt = (u & v).sum()
  1312. elif w is None:
  1313. ntt = (u * v).sum()
  1314. else:
  1315. w = _validate_weights(w)
  1316. ntt = (u * v * w).sum()
  1317. (nft, ntf) = _nbool_correspond_ft_tf(u, v, w=w)
  1318. denom = np.array(ntt + 2.0 * (ntf + nft))
  1319. if not denom.any():
  1320. raise ValueError('Sokal-Sneath dissimilarity is not defined for '
  1321. 'vectors that are entirely false.')
  1322. return float(2.0 * (ntf + nft)) / denom
  1323. _convert_to_double = partial(_convert_to_type, out_type=np.double)
  1324. _convert_to_bool = partial(_convert_to_type, out_type=bool)
  1325. # adding python-only wrappers to _distance_wrap module
  1326. _distance_wrap.pdist_correlation_double_wrap = _correlation_pdist_wrap
  1327. _distance_wrap.cdist_correlation_double_wrap = _correlation_cdist_wrap
  1328. @dataclasses.dataclass(frozen=True)
  1329. class CDistMetricWrapper:
  1330. metric_name: str
  1331. def __call__(self, XA, XB, *, out=None, **kwargs):
  1332. XA = np.ascontiguousarray(XA)
  1333. XB = np.ascontiguousarray(XB)
  1334. mA, n = XA.shape
  1335. mB, _ = XB.shape
  1336. metric_name = self.metric_name
  1337. metric_info = _METRICS[metric_name]
  1338. XA, XB, typ, kwargs = _validate_cdist_input(
  1339. XA, XB, mA, mB, n, metric_info, **kwargs)
  1340. w = kwargs.pop('w', None)
  1341. if w is not None:
  1342. metric = metric_info.dist_func
  1343. return _cdist_callable(
  1344. XA, XB, metric=metric, out=out, w=w, **kwargs)
  1345. dm = _prepare_out_argument(out, np.double, (mA, mB))
  1346. # get cdist wrapper
  1347. cdist_fn = getattr(_distance_wrap, f'cdist_{metric_name}_{typ}_wrap')
  1348. cdist_fn(XA, XB, dm, **kwargs)
  1349. return dm
  1350. @dataclasses.dataclass(frozen=True)
  1351. class CDistWeightedMetricWrapper:
  1352. metric_name: str
  1353. weighted_metric: str
  1354. def __call__(self, XA, XB, *, out=None, **kwargs):
  1355. XA = np.ascontiguousarray(XA)
  1356. XB = np.ascontiguousarray(XB)
  1357. mA, n = XA.shape
  1358. mB, _ = XB.shape
  1359. metric_name = self.metric_name
  1360. XA, XB, typ, kwargs = _validate_cdist_input(
  1361. XA, XB, mA, mB, n, _METRICS[metric_name], **kwargs)
  1362. dm = _prepare_out_argument(out, np.double, (mA, mB))
  1363. w = kwargs.pop('w', None)
  1364. if w is not None:
  1365. metric_name = self.weighted_metric
  1366. kwargs['w'] = w
  1367. # get cdist wrapper
  1368. cdist_fn = getattr(_distance_wrap, f'cdist_{metric_name}_{typ}_wrap')
  1369. cdist_fn(XA, XB, dm, **kwargs)
  1370. return dm
  1371. @dataclasses.dataclass(frozen=True)
  1372. class PDistMetricWrapper:
  1373. metric_name: str
  1374. def __call__(self, X, *, out=None, **kwargs):
  1375. X = np.ascontiguousarray(X)
  1376. m, n = X.shape
  1377. metric_name = self.metric_name
  1378. metric_info = _METRICS[metric_name]
  1379. X, typ, kwargs = _validate_pdist_input(
  1380. X, m, n, metric_info, **kwargs)
  1381. out_size = (m * (m - 1)) // 2
  1382. w = kwargs.pop('w', None)
  1383. if w is not None:
  1384. metric = metric_info.dist_func
  1385. return _pdist_callable(
  1386. X, metric=metric, out=out, w=w, **kwargs)
  1387. dm = _prepare_out_argument(out, np.double, (out_size,))
  1388. # get pdist wrapper
  1389. pdist_fn = getattr(_distance_wrap, f'pdist_{metric_name}_{typ}_wrap')
  1390. pdist_fn(X, dm, **kwargs)
  1391. return dm
  1392. @dataclasses.dataclass(frozen=True)
  1393. class PDistWeightedMetricWrapper:
  1394. metric_name: str
  1395. weighted_metric: str
  1396. def __call__(self, X, *, out=None, **kwargs):
  1397. X = np.ascontiguousarray(X)
  1398. m, n = X.shape
  1399. metric_name = self.metric_name
  1400. X, typ, kwargs = _validate_pdist_input(
  1401. X, m, n, _METRICS[metric_name], **kwargs)
  1402. out_size = (m * (m - 1)) // 2
  1403. dm = _prepare_out_argument(out, np.double, (out_size,))
  1404. w = kwargs.pop('w', None)
  1405. if w is not None:
  1406. metric_name = self.weighted_metric
  1407. kwargs['w'] = w
  1408. # get pdist wrapper
  1409. pdist_fn = getattr(_distance_wrap, f'pdist_{metric_name}_{typ}_wrap')
  1410. pdist_fn(X, dm, **kwargs)
  1411. return dm
  1412. @dataclasses.dataclass(frozen=True)
  1413. class MetricInfo:
  1414. # Name of python distance function
  1415. canonical_name: str
  1416. # All aliases, including canonical_name
  1417. aka: Set[str]
  1418. # unvectorized distance function
  1419. dist_func: Callable
  1420. # Optimized cdist function
  1421. cdist_func: Callable
  1422. # Optimized pdist function
  1423. pdist_func: Callable
  1424. # function that checks kwargs and computes default values:
  1425. # f(X, m, n, **kwargs)
  1426. validator: Optional[Callable] = None
  1427. # list of supported types:
  1428. # X (pdist) and XA (cdist) are used to choose the type. if there is no
  1429. # match the first type is used. Default double
  1430. types: List[str] = dataclasses.field(default_factory=lambda: ['double'])
  1431. # true if out array must be C-contiguous
  1432. requires_contiguous_out: bool = True
  1433. # Registry of implemented metrics:
  1434. _METRIC_INFOS = [
  1435. MetricInfo(
  1436. canonical_name='braycurtis',
  1437. aka={'braycurtis'},
  1438. dist_func=braycurtis,
  1439. cdist_func=_distance_pybind.cdist_braycurtis,
  1440. pdist_func=_distance_pybind.pdist_braycurtis,
  1441. ),
  1442. MetricInfo(
  1443. canonical_name='canberra',
  1444. aka={'canberra'},
  1445. dist_func=canberra,
  1446. cdist_func=_distance_pybind.cdist_canberra,
  1447. pdist_func=_distance_pybind.pdist_canberra,
  1448. ),
  1449. MetricInfo(
  1450. canonical_name='chebyshev',
  1451. aka={'chebychev', 'chebyshev', 'cheby', 'cheb', 'ch'},
  1452. dist_func=chebyshev,
  1453. cdist_func=_distance_pybind.cdist_chebyshev,
  1454. pdist_func=_distance_pybind.pdist_chebyshev,
  1455. ),
  1456. MetricInfo(
  1457. canonical_name='cityblock',
  1458. aka={'cityblock', 'cblock', 'cb', 'c'},
  1459. dist_func=cityblock,
  1460. cdist_func=_distance_pybind.cdist_cityblock,
  1461. pdist_func=_distance_pybind.pdist_cityblock,
  1462. ),
  1463. MetricInfo(
  1464. canonical_name='correlation',
  1465. aka={'correlation', 'co'},
  1466. dist_func=correlation,
  1467. cdist_func=CDistMetricWrapper('correlation'),
  1468. pdist_func=PDistMetricWrapper('correlation'),
  1469. ),
  1470. MetricInfo(
  1471. canonical_name='cosine',
  1472. aka={'cosine', 'cos'},
  1473. dist_func=cosine,
  1474. cdist_func=CDistMetricWrapper('cosine'),
  1475. pdist_func=PDistMetricWrapper('cosine'),
  1476. ),
  1477. MetricInfo(
  1478. canonical_name='dice',
  1479. aka={'dice'},
  1480. types=['bool'],
  1481. dist_func=dice,
  1482. cdist_func=CDistMetricWrapper('dice'),
  1483. pdist_func=PDistMetricWrapper('dice'),
  1484. ),
  1485. MetricInfo(
  1486. canonical_name='euclidean',
  1487. aka={'euclidean', 'euclid', 'eu', 'e'},
  1488. dist_func=euclidean,
  1489. cdist_func=_distance_pybind.cdist_euclidean,
  1490. pdist_func=_distance_pybind.pdist_euclidean,
  1491. ),
  1492. MetricInfo(
  1493. canonical_name='hamming',
  1494. aka={'matching', 'hamming', 'hamm', 'ha', 'h'},
  1495. types=['double', 'bool'],
  1496. validator=_validate_hamming_kwargs,
  1497. dist_func=hamming,
  1498. cdist_func=CDistWeightedMetricWrapper('hamming', 'hamming'),
  1499. pdist_func=PDistWeightedMetricWrapper('hamming', 'hamming'),
  1500. ),
  1501. MetricInfo(
  1502. canonical_name='jaccard',
  1503. aka={'jaccard', 'jacc', 'ja', 'j'},
  1504. types=['double', 'bool'],
  1505. dist_func=jaccard,
  1506. cdist_func=CDistMetricWrapper('jaccard'),
  1507. pdist_func=PDistMetricWrapper('jaccard'),
  1508. ),
  1509. MetricInfo(
  1510. canonical_name='jensenshannon',
  1511. aka={'jensenshannon', 'js'},
  1512. dist_func=jensenshannon,
  1513. cdist_func=CDistMetricWrapper('jensenshannon'),
  1514. pdist_func=PDistMetricWrapper('jensenshannon'),
  1515. ),
  1516. MetricInfo(
  1517. canonical_name='kulsinski',
  1518. aka={'kulsinski'},
  1519. types=['bool'],
  1520. dist_func=kulsinski,
  1521. cdist_func=CDistMetricWrapper('kulsinski'),
  1522. pdist_func=PDistMetricWrapper('kulsinski'),
  1523. ),
  1524. MetricInfo(
  1525. canonical_name='kulczynski1',
  1526. aka={'kulczynski1'},
  1527. types=['bool'],
  1528. dist_func=kulczynski1,
  1529. cdist_func=CDistMetricWrapper('kulczynski1'),
  1530. pdist_func=PDistMetricWrapper('kulczynski1'),
  1531. ),
  1532. MetricInfo(
  1533. canonical_name='mahalanobis',
  1534. aka={'mahalanobis', 'mahal', 'mah'},
  1535. validator=_validate_mahalanobis_kwargs,
  1536. dist_func=mahalanobis,
  1537. cdist_func=CDistMetricWrapper('mahalanobis'),
  1538. pdist_func=PDistMetricWrapper('mahalanobis'),
  1539. ),
  1540. MetricInfo(
  1541. canonical_name='minkowski',
  1542. aka={'minkowski', 'mi', 'm', 'pnorm'},
  1543. validator=_validate_minkowski_kwargs,
  1544. dist_func=minkowski,
  1545. cdist_func=_distance_pybind.cdist_minkowski,
  1546. pdist_func=_distance_pybind.pdist_minkowski,
  1547. ),
  1548. MetricInfo(
  1549. canonical_name='rogerstanimoto',
  1550. aka={'rogerstanimoto'},
  1551. types=['bool'],
  1552. dist_func=rogerstanimoto,
  1553. cdist_func=CDistMetricWrapper('rogerstanimoto'),
  1554. pdist_func=PDistMetricWrapper('rogerstanimoto'),
  1555. ),
  1556. MetricInfo(
  1557. canonical_name='russellrao',
  1558. aka={'russellrao'},
  1559. types=['bool'],
  1560. dist_func=russellrao,
  1561. cdist_func=CDistMetricWrapper('russellrao'),
  1562. pdist_func=PDistMetricWrapper('russellrao'),
  1563. ),
  1564. MetricInfo(
  1565. canonical_name='seuclidean',
  1566. aka={'seuclidean', 'se', 's'},
  1567. validator=_validate_seuclidean_kwargs,
  1568. dist_func=seuclidean,
  1569. cdist_func=CDistMetricWrapper('seuclidean'),
  1570. pdist_func=PDistMetricWrapper('seuclidean'),
  1571. ),
  1572. MetricInfo(
  1573. canonical_name='sokalmichener',
  1574. aka={'sokalmichener'},
  1575. types=['bool'],
  1576. dist_func=sokalmichener,
  1577. cdist_func=CDistMetricWrapper('sokalmichener'),
  1578. pdist_func=PDistMetricWrapper('sokalmichener'),
  1579. ),
  1580. MetricInfo(
  1581. canonical_name='sokalsneath',
  1582. aka={'sokalsneath'},
  1583. types=['bool'],
  1584. dist_func=sokalsneath,
  1585. cdist_func=CDistMetricWrapper('sokalsneath'),
  1586. pdist_func=PDistMetricWrapper('sokalsneath'),
  1587. ),
  1588. MetricInfo(
  1589. canonical_name='sqeuclidean',
  1590. aka={'sqeuclidean', 'sqe', 'sqeuclid'},
  1591. dist_func=sqeuclidean,
  1592. cdist_func=_distance_pybind.cdist_sqeuclidean,
  1593. pdist_func=_distance_pybind.pdist_sqeuclidean,
  1594. ),
  1595. MetricInfo(
  1596. canonical_name='yule',
  1597. aka={'yule'},
  1598. types=['bool'],
  1599. dist_func=yule,
  1600. cdist_func=CDistMetricWrapper('yule'),
  1601. pdist_func=PDistMetricWrapper('yule'),
  1602. ),
  1603. ]
  1604. _METRICS = {info.canonical_name: info for info in _METRIC_INFOS}
  1605. _METRIC_ALIAS = dict((alias, info)
  1606. for info in _METRIC_INFOS
  1607. for alias in info.aka)
  1608. _METRICS_NAMES = list(_METRICS.keys())
  1609. _TEST_METRICS = {'test_' + info.canonical_name: info for info in _METRIC_INFOS}
  1610. def pdist(X, metric='euclidean', *, out=None, **kwargs):
  1611. """
  1612. Pairwise distances between observations in n-dimensional space.
  1613. See Notes for common calling conventions.
  1614. Parameters
  1615. ----------
  1616. X : array_like
  1617. An m by n array of m original observations in an
  1618. n-dimensional space.
  1619. metric : str or function, optional
  1620. The distance metric to use. The distance function can
  1621. be 'braycurtis', 'canberra', 'chebyshev', 'cityblock',
  1622. 'correlation', 'cosine', 'dice', 'euclidean', 'hamming',
  1623. 'jaccard', 'jensenshannon', 'kulczynski1',
  1624. 'mahalanobis', 'matching', 'minkowski', 'rogerstanimoto',
  1625. 'russellrao', 'seuclidean', 'sokalmichener', 'sokalsneath',
  1626. 'sqeuclidean', 'yule'.
  1627. **kwargs : dict, optional
  1628. Extra arguments to `metric`: refer to each metric documentation for a
  1629. list of all possible arguments.
  1630. Some possible arguments:
  1631. p : scalar
  1632. The p-norm to apply for Minkowski, weighted and unweighted.
  1633. Default: 2.
  1634. w : ndarray
  1635. The weight vector for metrics that support weights (e.g., Minkowski).
  1636. V : ndarray
  1637. The variance vector for standardized Euclidean.
  1638. Default: var(X, axis=0, ddof=1)
  1639. VI : ndarray
  1640. The inverse of the covariance matrix for Mahalanobis.
  1641. Default: inv(cov(X.T)).T
  1642. out : ndarray.
  1643. The output array
  1644. If not None, condensed distance matrix Y is stored in this array.
  1645. Returns
  1646. -------
  1647. Y : ndarray
  1648. Returns a condensed distance matrix Y. For each :math:`i` and :math:`j`
  1649. (where :math:`i<j<m`),where m is the number of original observations.
  1650. The metric ``dist(u=X[i], v=X[j])`` is computed and stored in entry ``m
  1651. * i + j - ((i + 2) * (i + 1)) // 2``.
  1652. See Also
  1653. --------
  1654. squareform : converts between condensed distance matrices and
  1655. square distance matrices.
  1656. Notes
  1657. -----
  1658. See ``squareform`` for information on how to calculate the index of
  1659. this entry or to convert the condensed distance matrix to a
  1660. redundant square matrix.
  1661. The following are common calling conventions.
  1662. 1. ``Y = pdist(X, 'euclidean')``
  1663. Computes the distance between m points using Euclidean distance
  1664. (2-norm) as the distance metric between the points. The points
  1665. are arranged as m n-dimensional row vectors in the matrix X.
  1666. 2. ``Y = pdist(X, 'minkowski', p=2.)``
  1667. Computes the distances using the Minkowski distance
  1668. :math:`\\|u-v\\|_p` (:math:`p`-norm) where :math:`p > 0` (note
  1669. that this is only a quasi-metric if :math:`0 < p < 1`).
  1670. 3. ``Y = pdist(X, 'cityblock')``
  1671. Computes the city block or Manhattan distance between the
  1672. points.
  1673. 4. ``Y = pdist(X, 'seuclidean', V=None)``
  1674. Computes the standardized Euclidean distance. The standardized
  1675. Euclidean distance between two n-vectors ``u`` and ``v`` is
  1676. .. math::
  1677. \\sqrt{\\sum {(u_i-v_i)^2 / V[x_i]}}
  1678. V is the variance vector; V[i] is the variance computed over all
  1679. the i'th components of the points. If not passed, it is
  1680. automatically computed.
  1681. 5. ``Y = pdist(X, 'sqeuclidean')``
  1682. Computes the squared Euclidean distance :math:`\\|u-v\\|_2^2` between
  1683. the vectors.
  1684. 6. ``Y = pdist(X, 'cosine')``
  1685. Computes the cosine distance between vectors u and v,
  1686. .. math::
  1687. 1 - \\frac{u \\cdot v}
  1688. {{\\|u\\|}_2 {\\|v\\|}_2}
  1689. where :math:`\\|*\\|_2` is the 2-norm of its argument ``*``, and
  1690. :math:`u \\cdot v` is the dot product of ``u`` and ``v``.
  1691. 7. ``Y = pdist(X, 'correlation')``
  1692. Computes the correlation distance between vectors u and v. This is
  1693. .. math::
  1694. 1 - \\frac{(u - \\bar{u}) \\cdot (v - \\bar{v})}
  1695. {{\\|(u - \\bar{u})\\|}_2 {\\|(v - \\bar{v})\\|}_2}
  1696. where :math:`\\bar{v}` is the mean of the elements of vector v,
  1697. and :math:`x \\cdot y` is the dot product of :math:`x` and :math:`y`.
  1698. 8. ``Y = pdist(X, 'hamming')``
  1699. Computes the normalized Hamming distance, or the proportion of
  1700. those vector elements between two n-vectors ``u`` and ``v``
  1701. which disagree. To save memory, the matrix ``X`` can be of type
  1702. boolean.
  1703. 9. ``Y = pdist(X, 'jaccard')``
  1704. Computes the Jaccard distance between the points. Given two
  1705. vectors, ``u`` and ``v``, the Jaccard distance is the
  1706. proportion of those elements ``u[i]`` and ``v[i]`` that
  1707. disagree.
  1708. 10. ``Y = pdist(X, 'jensenshannon')``
  1709. Computes the Jensen-Shannon distance between two probability arrays.
  1710. Given two probability vectors, :math:`p` and :math:`q`, the
  1711. Jensen-Shannon distance is
  1712. .. math::
  1713. \\sqrt{\\frac{D(p \\parallel m) + D(q \\parallel m)}{2}}
  1714. where :math:`m` is the pointwise mean of :math:`p` and :math:`q`
  1715. and :math:`D` is the Kullback-Leibler divergence.
  1716. 11. ``Y = pdist(X, 'chebyshev')``
  1717. Computes the Chebyshev distance between the points. The
  1718. Chebyshev distance between two n-vectors ``u`` and ``v`` is the
  1719. maximum norm-1 distance between their respective elements. More
  1720. precisely, the distance is given by
  1721. .. math::
  1722. d(u,v) = \\max_i {|u_i-v_i|}
  1723. 12. ``Y = pdist(X, 'canberra')``
  1724. Computes the Canberra distance between the points. The
  1725. Canberra distance between two points ``u`` and ``v`` is
  1726. .. math::
  1727. d(u,v) = \\sum_i \\frac{|u_i-v_i|}
  1728. {|u_i|+|v_i|}
  1729. 13. ``Y = pdist(X, 'braycurtis')``
  1730. Computes the Bray-Curtis distance between the points. The
  1731. Bray-Curtis distance between two points ``u`` and ``v`` is
  1732. .. math::
  1733. d(u,v) = \\frac{\\sum_i {|u_i-v_i|}}
  1734. {\\sum_i {|u_i+v_i|}}
  1735. 14. ``Y = pdist(X, 'mahalanobis', VI=None)``
  1736. Computes the Mahalanobis distance between the points. The
  1737. Mahalanobis distance between two points ``u`` and ``v`` is
  1738. :math:`\\sqrt{(u-v)(1/V)(u-v)^T}` where :math:`(1/V)` (the ``VI``
  1739. variable) is the inverse covariance. If ``VI`` is not None,
  1740. ``VI`` will be used as the inverse covariance matrix.
  1741. 15. ``Y = pdist(X, 'yule')``
  1742. Computes the Yule distance between each pair of boolean
  1743. vectors. (see yule function documentation)
  1744. 16. ``Y = pdist(X, 'matching')``
  1745. Synonym for 'hamming'.
  1746. 17. ``Y = pdist(X, 'dice')``
  1747. Computes the Dice distance between each pair of boolean
  1748. vectors. (see dice function documentation)
  1749. 18. ``Y = pdist(X, 'kulczynski1')``
  1750. Computes the kulczynski1 distance between each pair of
  1751. boolean vectors. (see kulczynski1 function documentation)
  1752. 19. ``Y = pdist(X, 'rogerstanimoto')``
  1753. Computes the Rogers-Tanimoto distance between each pair of
  1754. boolean vectors. (see rogerstanimoto function documentation)
  1755. 20. ``Y = pdist(X, 'russellrao')``
  1756. Computes the Russell-Rao distance between each pair of
  1757. boolean vectors. (see russellrao function documentation)
  1758. 21. ``Y = pdist(X, 'sokalmichener')``
  1759. Computes the Sokal-Michener distance between each pair of
  1760. boolean vectors. (see sokalmichener function documentation)
  1761. 22. ``Y = pdist(X, 'sokalsneath')``
  1762. Computes the Sokal-Sneath distance between each pair of
  1763. boolean vectors. (see sokalsneath function documentation)
  1764. 23. ``Y = pdist(X, 'kulczynski1')``
  1765. Computes the Kulczynski 1 distance between each pair of
  1766. boolean vectors. (see kulczynski1 function documentation)
  1767. 24. ``Y = pdist(X, f)``
  1768. Computes the distance between all pairs of vectors in X
  1769. using the user supplied 2-arity function f. For example,
  1770. Euclidean distance between the vectors could be computed
  1771. as follows::
  1772. dm = pdist(X, lambda u, v: np.sqrt(((u-v)**2).sum()))
  1773. Note that you should avoid passing a reference to one of
  1774. the distance functions defined in this library. For example,::
  1775. dm = pdist(X, sokalsneath)
  1776. would calculate the pair-wise distances between the vectors in
  1777. X using the Python function sokalsneath. This would result in
  1778. sokalsneath being called :math:`{n \\choose 2}` times, which
  1779. is inefficient. Instead, the optimized C version is more
  1780. efficient, and we call it using the following syntax.::
  1781. dm = pdist(X, 'sokalsneath')
  1782. """
  1783. # You can also call this as:
  1784. # Y = pdist(X, 'test_abc')
  1785. # where 'abc' is the metric being tested. This computes the distance
  1786. # between all pairs of vectors in X using the distance metric 'abc' but
  1787. # with a more succinct, verifiable, but less efficient implementation.
  1788. X = _asarray_validated(X, sparse_ok=False, objects_ok=True, mask_ok=True,
  1789. check_finite=False)
  1790. s = X.shape
  1791. if len(s) != 2:
  1792. raise ValueError('A 2-dimensional array must be passed.')
  1793. m, n = s
  1794. if callable(metric):
  1795. mstr = getattr(metric, '__name__', 'UnknownCustomMetric')
  1796. metric_info = _METRIC_ALIAS.get(mstr, None)
  1797. if metric_info is not None:
  1798. X, typ, kwargs = _validate_pdist_input(
  1799. X, m, n, metric_info, **kwargs)
  1800. return _pdist_callable(X, metric=metric, out=out, **kwargs)
  1801. elif isinstance(metric, str):
  1802. mstr = metric.lower()
  1803. metric_info = _METRIC_ALIAS.get(mstr, None)
  1804. if metric_info is not None:
  1805. pdist_fn = metric_info.pdist_func
  1806. return pdist_fn(X, out=out, **kwargs)
  1807. elif mstr.startswith("test_"):
  1808. metric_info = _TEST_METRICS.get(mstr, None)
  1809. if metric_info is None:
  1810. raise ValueError(f'Unknown "Test" Distance Metric: {mstr[5:]}')
  1811. X, typ, kwargs = _validate_pdist_input(
  1812. X, m, n, metric_info, **kwargs)
  1813. return _pdist_callable(
  1814. X, metric=metric_info.dist_func, out=out, **kwargs)
  1815. else:
  1816. raise ValueError('Unknown Distance Metric: %s' % mstr)
  1817. else:
  1818. raise TypeError('2nd argument metric must be a string identifier '
  1819. 'or a function.')
  1820. def squareform(X, force="no", checks=True):
  1821. """
  1822. Convert a vector-form distance vector to a square-form distance
  1823. matrix, and vice-versa.
  1824. Parameters
  1825. ----------
  1826. X : array_like
  1827. Either a condensed or redundant distance matrix.
  1828. force : str, optional
  1829. As with MATLAB(TM), if force is equal to ``'tovector'`` or
  1830. ``'tomatrix'``, the input will be treated as a distance matrix or
  1831. distance vector respectively.
  1832. checks : bool, optional
  1833. If set to False, no checks will be made for matrix
  1834. symmetry nor zero diagonals. This is useful if it is known that
  1835. ``X - X.T1`` is small and ``diag(X)`` is close to zero.
  1836. These values are ignored any way so they do not disrupt the
  1837. squareform transformation.
  1838. Returns
  1839. -------
  1840. Y : ndarray
  1841. If a condensed distance matrix is passed, a redundant one is
  1842. returned, or if a redundant one is passed, a condensed distance
  1843. matrix is returned.
  1844. Notes
  1845. -----
  1846. 1. ``v = squareform(X)``
  1847. Given a square n-by-n symmetric distance matrix ``X``,
  1848. ``v = squareform(X)`` returns a ``n * (n-1) / 2``
  1849. (i.e. binomial coefficient n choose 2) sized vector `v`
  1850. where :math:`v[{n \\choose 2} - {n-i \\choose 2} + (j-i-1)]`
  1851. is the distance between distinct points ``i`` and ``j``.
  1852. If ``X`` is non-square or asymmetric, an error is raised.
  1853. 2. ``X = squareform(v)``
  1854. Given a ``n * (n-1) / 2`` sized vector ``v``
  1855. for some integer ``n >= 1`` encoding distances as described,
  1856. ``X = squareform(v)`` returns a n-by-n distance matrix ``X``.
  1857. The ``X[i, j]`` and ``X[j, i]`` values are set to
  1858. :math:`v[{n \\choose 2} - {n-i \\choose 2} + (j-i-1)]`
  1859. and all diagonal elements are zero.
  1860. In SciPy 0.19.0, ``squareform`` stopped casting all input types to
  1861. float64, and started returning arrays of the same dtype as the input.
  1862. """
  1863. X = np.ascontiguousarray(X)
  1864. s = X.shape
  1865. if force.lower() == 'tomatrix':
  1866. if len(s) != 1:
  1867. raise ValueError("Forcing 'tomatrix' but input X is not a "
  1868. "distance vector.")
  1869. elif force.lower() == 'tovector':
  1870. if len(s) != 2:
  1871. raise ValueError("Forcing 'tovector' but input X is not a "
  1872. "distance matrix.")
  1873. # X = squareform(v)
  1874. if len(s) == 1:
  1875. if s[0] == 0:
  1876. return np.zeros((1, 1), dtype=X.dtype)
  1877. # Grab the closest value to the square root of the number
  1878. # of elements times 2 to see if the number of elements
  1879. # is indeed a binomial coefficient.
  1880. d = int(np.ceil(np.sqrt(s[0] * 2)))
  1881. # Check that v is of valid dimensions.
  1882. if d * (d - 1) != s[0] * 2:
  1883. raise ValueError('Incompatible vector size. It must be a binomial '
  1884. 'coefficient n choose 2 for some integer n >= 2.')
  1885. # Allocate memory for the distance matrix.
  1886. M = np.zeros((d, d), dtype=X.dtype)
  1887. # Since the C code does not support striding using strides.
  1888. # The dimensions are used instead.
  1889. X = _copy_array_if_base_present(X)
  1890. # Fill in the values of the distance matrix.
  1891. _distance_wrap.to_squareform_from_vector_wrap(M, X)
  1892. # Return the distance matrix.
  1893. return M
  1894. elif len(s) == 2:
  1895. if s[0] != s[1]:
  1896. raise ValueError('The matrix argument must be square.')
  1897. if checks:
  1898. is_valid_dm(X, throw=True, name='X')
  1899. # One-side of the dimensions is set here.
  1900. d = s[0]
  1901. if d <= 1:
  1902. return np.array([], dtype=X.dtype)
  1903. # Create a vector.
  1904. v = np.zeros((d * (d - 1)) // 2, dtype=X.dtype)
  1905. # Since the C code does not support striding using strides.
  1906. # The dimensions are used instead.
  1907. X = _copy_array_if_base_present(X)
  1908. # Convert the vector to squareform.
  1909. _distance_wrap.to_vector_from_squareform_wrap(X, v)
  1910. return v
  1911. else:
  1912. raise ValueError(('The first argument must be one or two dimensional '
  1913. 'array. A %d-dimensional array is not '
  1914. 'permitted') % len(s))
  1915. def is_valid_dm(D, tol=0.0, throw=False, name="D", warning=False):
  1916. """
  1917. Return True if input array is a valid distance matrix.
  1918. Distance matrices must be 2-dimensional numpy arrays.
  1919. They must have a zero-diagonal, and they must be symmetric.
  1920. Parameters
  1921. ----------
  1922. D : array_like
  1923. The candidate object to test for validity.
  1924. tol : float, optional
  1925. The distance matrix should be symmetric. `tol` is the maximum
  1926. difference between entries ``ij`` and ``ji`` for the distance
  1927. metric to be considered symmetric.
  1928. throw : bool, optional
  1929. An exception is thrown if the distance matrix passed is not valid.
  1930. name : str, optional
  1931. The name of the variable to checked. This is useful if
  1932. throw is set to True so the offending variable can be identified
  1933. in the exception message when an exception is thrown.
  1934. warning : bool, optional
  1935. Instead of throwing an exception, a warning message is
  1936. raised.
  1937. Returns
  1938. -------
  1939. valid : bool
  1940. True if the variable `D` passed is a valid distance matrix.
  1941. Notes
  1942. -----
  1943. Small numerical differences in `D` and `D.T` and non-zeroness of
  1944. the diagonal are ignored if they are within the tolerance specified
  1945. by `tol`.
  1946. """
  1947. D = np.asarray(D, order='c')
  1948. valid = True
  1949. try:
  1950. s = D.shape
  1951. if len(D.shape) != 2:
  1952. if name:
  1953. raise ValueError(('Distance matrix \'%s\' must have shape=2 '
  1954. '(i.e. be two-dimensional).') % name)
  1955. else:
  1956. raise ValueError('Distance matrix must have shape=2 (i.e. '
  1957. 'be two-dimensional).')
  1958. if tol == 0.0:
  1959. if not (D == D.T).all():
  1960. if name:
  1961. raise ValueError(('Distance matrix \'%s\' must be '
  1962. 'symmetric.') % name)
  1963. else:
  1964. raise ValueError('Distance matrix must be symmetric.')
  1965. if not (D[range(0, s[0]), range(0, s[0])] == 0).all():
  1966. if name:
  1967. raise ValueError(('Distance matrix \'%s\' diagonal must '
  1968. 'be zero.') % name)
  1969. else:
  1970. raise ValueError('Distance matrix diagonal must be zero.')
  1971. else:
  1972. if not (D - D.T <= tol).all():
  1973. if name:
  1974. raise ValueError(('Distance matrix \'%s\' must be '
  1975. 'symmetric within tolerance %5.5f.')
  1976. % (name, tol))
  1977. else:
  1978. raise ValueError('Distance matrix must be symmetric within'
  1979. ' tolerance %5.5f.' % tol)
  1980. if not (D[range(0, s[0]), range(0, s[0])] <= tol).all():
  1981. if name:
  1982. raise ValueError(('Distance matrix \'%s\' diagonal must be'
  1983. ' close to zero within tolerance %5.5f.')
  1984. % (name, tol))
  1985. else:
  1986. raise ValueError(('Distance matrix \'%s\' diagonal must be'
  1987. ' close to zero within tolerance %5.5f.')
  1988. % tol)
  1989. except Exception as e:
  1990. if throw:
  1991. raise
  1992. if warning:
  1993. warnings.warn(str(e))
  1994. valid = False
  1995. return valid
  1996. def is_valid_y(y, warning=False, throw=False, name=None):
  1997. """
  1998. Return True if the input array is a valid condensed distance matrix.
  1999. Condensed distance matrices must be 1-dimensional numpy arrays.
  2000. Their length must be a binomial coefficient :math:`{n \\choose 2}`
  2001. for some positive integer n.
  2002. Parameters
  2003. ----------
  2004. y : array_like
  2005. The condensed distance matrix.
  2006. warning : bool, optional
  2007. Invokes a warning if the variable passed is not a valid
  2008. condensed distance matrix. The warning message explains why
  2009. the distance matrix is not valid. `name` is used when
  2010. referencing the offending variable.
  2011. throw : bool, optional
  2012. Throws an exception if the variable passed is not a valid
  2013. condensed distance matrix.
  2014. name : bool, optional
  2015. Used when referencing the offending variable in the
  2016. warning or exception message.
  2017. """
  2018. y = np.asarray(y, order='c')
  2019. valid = True
  2020. try:
  2021. if len(y.shape) != 1:
  2022. if name:
  2023. raise ValueError(('Condensed distance matrix \'%s\' must '
  2024. 'have shape=1 (i.e. be one-dimensional).')
  2025. % name)
  2026. else:
  2027. raise ValueError('Condensed distance matrix must have shape=1 '
  2028. '(i.e. be one-dimensional).')
  2029. n = y.shape[0]
  2030. d = int(np.ceil(np.sqrt(n * 2)))
  2031. if (d * (d - 1) / 2) != n:
  2032. if name:
  2033. raise ValueError(('Length n of condensed distance matrix '
  2034. '\'%s\' must be a binomial coefficient, i.e.'
  2035. 'there must be a k such that '
  2036. '(k \\choose 2)=n)!') % name)
  2037. else:
  2038. raise ValueError('Length n of condensed distance matrix must '
  2039. 'be a binomial coefficient, i.e. there must '
  2040. 'be a k such that (k \\choose 2)=n)!')
  2041. except Exception as e:
  2042. if throw:
  2043. raise
  2044. if warning:
  2045. warnings.warn(str(e))
  2046. valid = False
  2047. return valid
  2048. def num_obs_dm(d):
  2049. """
  2050. Return the number of original observations that correspond to a
  2051. square, redundant distance matrix.
  2052. Parameters
  2053. ----------
  2054. d : array_like
  2055. The target distance matrix.
  2056. Returns
  2057. -------
  2058. num_obs_dm : int
  2059. The number of observations in the redundant distance matrix.
  2060. """
  2061. d = np.asarray(d, order='c')
  2062. is_valid_dm(d, tol=np.inf, throw=True, name='d')
  2063. return d.shape[0]
  2064. def num_obs_y(Y):
  2065. """
  2066. Return the number of original observations that correspond to a
  2067. condensed distance matrix.
  2068. Parameters
  2069. ----------
  2070. Y : array_like
  2071. Condensed distance matrix.
  2072. Returns
  2073. -------
  2074. n : int
  2075. The number of observations in the condensed distance matrix `Y`.
  2076. """
  2077. Y = np.asarray(Y, order='c')
  2078. is_valid_y(Y, throw=True, name='Y')
  2079. k = Y.shape[0]
  2080. if k == 0:
  2081. raise ValueError("The number of observations cannot be determined on "
  2082. "an empty distance matrix.")
  2083. d = int(np.ceil(np.sqrt(k * 2)))
  2084. if (d * (d - 1) / 2) != k:
  2085. raise ValueError("Invalid condensed distance matrix passed. Must be "
  2086. "some k where k=(n choose 2) for some n >= 2.")
  2087. return d
  2088. def _prepare_out_argument(out, dtype, expected_shape):
  2089. if out is None:
  2090. return np.empty(expected_shape, dtype=dtype)
  2091. if out.shape != expected_shape:
  2092. raise ValueError("Output array has incorrect shape.")
  2093. if not out.flags.c_contiguous:
  2094. raise ValueError("Output array must be C-contiguous.")
  2095. if out.dtype != np.double:
  2096. raise ValueError("Output array must be double type.")
  2097. return out
  2098. def _pdist_callable(X, *, out, metric, **kwargs):
  2099. n = X.shape[0]
  2100. out_size = (n * (n - 1)) // 2
  2101. dm = _prepare_out_argument(out, np.double, (out_size,))
  2102. k = 0
  2103. for i in range(X.shape[0] - 1):
  2104. for j in range(i + 1, X.shape[0]):
  2105. dm[k] = metric(X[i], X[j], **kwargs)
  2106. k += 1
  2107. return dm
  2108. def _cdist_callable(XA, XB, *, out, metric, **kwargs):
  2109. mA = XA.shape[0]
  2110. mB = XB.shape[0]
  2111. dm = _prepare_out_argument(out, np.double, (mA, mB))
  2112. for i in range(mA):
  2113. for j in range(mB):
  2114. dm[i, j] = metric(XA[i], XB[j], **kwargs)
  2115. return dm
  2116. def cdist(XA, XB, metric='euclidean', *, out=None, **kwargs):
  2117. """
  2118. Compute distance between each pair of the two collections of inputs.
  2119. See Notes for common calling conventions.
  2120. Parameters
  2121. ----------
  2122. XA : array_like
  2123. An :math:`m_A` by :math:`n` array of :math:`m_A`
  2124. original observations in an :math:`n`-dimensional space.
  2125. Inputs are converted to float type.
  2126. XB : array_like
  2127. An :math:`m_B` by :math:`n` array of :math:`m_B`
  2128. original observations in an :math:`n`-dimensional space.
  2129. Inputs are converted to float type.
  2130. metric : str or callable, optional
  2131. The distance metric to use. If a string, the distance function can be
  2132. 'braycurtis', 'canberra', 'chebyshev', 'cityblock', 'correlation',
  2133. 'cosine', 'dice', 'euclidean', 'hamming', 'jaccard', 'jensenshannon',
  2134. 'kulczynski1', 'mahalanobis', 'matching', 'minkowski',
  2135. 'rogerstanimoto', 'russellrao', 'seuclidean', 'sokalmichener',
  2136. 'sokalsneath', 'sqeuclidean', 'yule'.
  2137. **kwargs : dict, optional
  2138. Extra arguments to `metric`: refer to each metric documentation for a
  2139. list of all possible arguments.
  2140. Some possible arguments:
  2141. p : scalar
  2142. The p-norm to apply for Minkowski, weighted and unweighted.
  2143. Default: 2.
  2144. w : array_like
  2145. The weight vector for metrics that support weights (e.g., Minkowski).
  2146. V : array_like
  2147. The variance vector for standardized Euclidean.
  2148. Default: var(vstack([XA, XB]), axis=0, ddof=1)
  2149. VI : array_like
  2150. The inverse of the covariance matrix for Mahalanobis.
  2151. Default: inv(cov(vstack([XA, XB].T))).T
  2152. out : ndarray
  2153. The output array
  2154. If not None, the distance matrix Y is stored in this array.
  2155. Returns
  2156. -------
  2157. Y : ndarray
  2158. A :math:`m_A` by :math:`m_B` distance matrix is returned.
  2159. For each :math:`i` and :math:`j`, the metric
  2160. ``dist(u=XA[i], v=XB[j])`` is computed and stored in the
  2161. :math:`ij` th entry.
  2162. Raises
  2163. ------
  2164. ValueError
  2165. An exception is thrown if `XA` and `XB` do not have
  2166. the same number of columns.
  2167. Notes
  2168. -----
  2169. The following are common calling conventions:
  2170. 1. ``Y = cdist(XA, XB, 'euclidean')``
  2171. Computes the distance between :math:`m` points using
  2172. Euclidean distance (2-norm) as the distance metric between the
  2173. points. The points are arranged as :math:`m`
  2174. :math:`n`-dimensional row vectors in the matrix X.
  2175. 2. ``Y = cdist(XA, XB, 'minkowski', p=2.)``
  2176. Computes the distances using the Minkowski distance
  2177. :math:`\\|u-v\\|_p` (:math:`p`-norm) where :math:`p > 0` (note
  2178. that this is only a quasi-metric if :math:`0 < p < 1`).
  2179. 3. ``Y = cdist(XA, XB, 'cityblock')``
  2180. Computes the city block or Manhattan distance between the
  2181. points.
  2182. 4. ``Y = cdist(XA, XB, 'seuclidean', V=None)``
  2183. Computes the standardized Euclidean distance. The standardized
  2184. Euclidean distance between two n-vectors ``u`` and ``v`` is
  2185. .. math::
  2186. \\sqrt{\\sum {(u_i-v_i)^2 / V[x_i]}}.
  2187. V is the variance vector; V[i] is the variance computed over all
  2188. the i'th components of the points. If not passed, it is
  2189. automatically computed.
  2190. 5. ``Y = cdist(XA, XB, 'sqeuclidean')``
  2191. Computes the squared Euclidean distance :math:`\\|u-v\\|_2^2` between
  2192. the vectors.
  2193. 6. ``Y = cdist(XA, XB, 'cosine')``
  2194. Computes the cosine distance between vectors u and v,
  2195. .. math::
  2196. 1 - \\frac{u \\cdot v}
  2197. {{\\|u\\|}_2 {\\|v\\|}_2}
  2198. where :math:`\\|*\\|_2` is the 2-norm of its argument ``*``, and
  2199. :math:`u \\cdot v` is the dot product of :math:`u` and :math:`v`.
  2200. 7. ``Y = cdist(XA, XB, 'correlation')``
  2201. Computes the correlation distance between vectors u and v. This is
  2202. .. math::
  2203. 1 - \\frac{(u - \\bar{u}) \\cdot (v - \\bar{v})}
  2204. {{\\|(u - \\bar{u})\\|}_2 {\\|(v - \\bar{v})\\|}_2}
  2205. where :math:`\\bar{v}` is the mean of the elements of vector v,
  2206. and :math:`x \\cdot y` is the dot product of :math:`x` and :math:`y`.
  2207. 8. ``Y = cdist(XA, XB, 'hamming')``
  2208. Computes the normalized Hamming distance, or the proportion of
  2209. those vector elements between two n-vectors ``u`` and ``v``
  2210. which disagree. To save memory, the matrix ``X`` can be of type
  2211. boolean.
  2212. 9. ``Y = cdist(XA, XB, 'jaccard')``
  2213. Computes the Jaccard distance between the points. Given two
  2214. vectors, ``u`` and ``v``, the Jaccard distance is the
  2215. proportion of those elements ``u[i]`` and ``v[i]`` that
  2216. disagree where at least one of them is non-zero.
  2217. 10. ``Y = cdist(XA, XB, 'jensenshannon')``
  2218. Computes the Jensen-Shannon distance between two probability arrays.
  2219. Given two probability vectors, :math:`p` and :math:`q`, the
  2220. Jensen-Shannon distance is
  2221. .. math::
  2222. \\sqrt{\\frac{D(p \\parallel m) + D(q \\parallel m)}{2}}
  2223. where :math:`m` is the pointwise mean of :math:`p` and :math:`q`
  2224. and :math:`D` is the Kullback-Leibler divergence.
  2225. 11. ``Y = cdist(XA, XB, 'chebyshev')``
  2226. Computes the Chebyshev distance between the points. The
  2227. Chebyshev distance between two n-vectors ``u`` and ``v`` is the
  2228. maximum norm-1 distance between their respective elements. More
  2229. precisely, the distance is given by
  2230. .. math::
  2231. d(u,v) = \\max_i {|u_i-v_i|}.
  2232. 12. ``Y = cdist(XA, XB, 'canberra')``
  2233. Computes the Canberra distance between the points. The
  2234. Canberra distance between two points ``u`` and ``v`` is
  2235. .. math::
  2236. d(u,v) = \\sum_i \\frac{|u_i-v_i|}
  2237. {|u_i|+|v_i|}.
  2238. 13. ``Y = cdist(XA, XB, 'braycurtis')``
  2239. Computes the Bray-Curtis distance between the points. The
  2240. Bray-Curtis distance between two points ``u`` and ``v`` is
  2241. .. math::
  2242. d(u,v) = \\frac{\\sum_i (|u_i-v_i|)}
  2243. {\\sum_i (|u_i+v_i|)}
  2244. 14. ``Y = cdist(XA, XB, 'mahalanobis', VI=None)``
  2245. Computes the Mahalanobis distance between the points. The
  2246. Mahalanobis distance between two points ``u`` and ``v`` is
  2247. :math:`\\sqrt{(u-v)(1/V)(u-v)^T}` where :math:`(1/V)` (the ``VI``
  2248. variable) is the inverse covariance. If ``VI`` is not None,
  2249. ``VI`` will be used as the inverse covariance matrix.
  2250. 15. ``Y = cdist(XA, XB, 'yule')``
  2251. Computes the Yule distance between the boolean
  2252. vectors. (see `yule` function documentation)
  2253. 16. ``Y = cdist(XA, XB, 'matching')``
  2254. Synonym for 'hamming'.
  2255. 17. ``Y = cdist(XA, XB, 'dice')``
  2256. Computes the Dice distance between the boolean vectors. (see
  2257. `dice` function documentation)
  2258. 18. ``Y = cdist(XA, XB, 'kulczynski1')``
  2259. Computes the kulczynski distance between the boolean
  2260. vectors. (see `kulczynski1` function documentation)
  2261. 19. ``Y = cdist(XA, XB, 'rogerstanimoto')``
  2262. Computes the Rogers-Tanimoto distance between the boolean
  2263. vectors. (see `rogerstanimoto` function documentation)
  2264. 20. ``Y = cdist(XA, XB, 'russellrao')``
  2265. Computes the Russell-Rao distance between the boolean
  2266. vectors. (see `russellrao` function documentation)
  2267. 21. ``Y = cdist(XA, XB, 'sokalmichener')``
  2268. Computes the Sokal-Michener distance between the boolean
  2269. vectors. (see `sokalmichener` function documentation)
  2270. 22. ``Y = cdist(XA, XB, 'sokalsneath')``
  2271. Computes the Sokal-Sneath distance between the vectors. (see
  2272. `sokalsneath` function documentation)
  2273. 23. ``Y = cdist(XA, XB, f)``
  2274. Computes the distance between all pairs of vectors in X
  2275. using the user supplied 2-arity function f. For example,
  2276. Euclidean distance between the vectors could be computed
  2277. as follows::
  2278. dm = cdist(XA, XB, lambda u, v: np.sqrt(((u-v)**2).sum()))
  2279. Note that you should avoid passing a reference to one of
  2280. the distance functions defined in this library. For example,::
  2281. dm = cdist(XA, XB, sokalsneath)
  2282. would calculate the pair-wise distances between the vectors in
  2283. X using the Python function `sokalsneath`. This would result in
  2284. sokalsneath being called :math:`{n \\choose 2}` times, which
  2285. is inefficient. Instead, the optimized C version is more
  2286. efficient, and we call it using the following syntax::
  2287. dm = cdist(XA, XB, 'sokalsneath')
  2288. Examples
  2289. --------
  2290. Find the Euclidean distances between four 2-D coordinates:
  2291. >>> from scipy.spatial import distance
  2292. >>> import numpy as np
  2293. >>> coords = [(35.0456, -85.2672),
  2294. ... (35.1174, -89.9711),
  2295. ... (35.9728, -83.9422),
  2296. ... (36.1667, -86.7833)]
  2297. >>> distance.cdist(coords, coords, 'euclidean')
  2298. array([[ 0. , 4.7044, 1.6172, 1.8856],
  2299. [ 4.7044, 0. , 6.0893, 3.3561],
  2300. [ 1.6172, 6.0893, 0. , 2.8477],
  2301. [ 1.8856, 3.3561, 2.8477, 0. ]])
  2302. Find the Manhattan distance from a 3-D point to the corners of the unit
  2303. cube:
  2304. >>> a = np.array([[0, 0, 0],
  2305. ... [0, 0, 1],
  2306. ... [0, 1, 0],
  2307. ... [0, 1, 1],
  2308. ... [1, 0, 0],
  2309. ... [1, 0, 1],
  2310. ... [1, 1, 0],
  2311. ... [1, 1, 1]])
  2312. >>> b = np.array([[ 0.1, 0.2, 0.4]])
  2313. >>> distance.cdist(a, b, 'cityblock')
  2314. array([[ 0.7],
  2315. [ 0.9],
  2316. [ 1.3],
  2317. [ 1.5],
  2318. [ 1.5],
  2319. [ 1.7],
  2320. [ 2.1],
  2321. [ 2.3]])
  2322. """
  2323. # You can also call this as:
  2324. # Y = cdist(XA, XB, 'test_abc')
  2325. # where 'abc' is the metric being tested. This computes the distance
  2326. # between all pairs of vectors in XA and XB using the distance metric 'abc'
  2327. # but with a more succinct, verifiable, but less efficient implementation.
  2328. XA = np.asarray(XA)
  2329. XB = np.asarray(XB)
  2330. s = XA.shape
  2331. sB = XB.shape
  2332. if len(s) != 2:
  2333. raise ValueError('XA must be a 2-dimensional array.')
  2334. if len(sB) != 2:
  2335. raise ValueError('XB must be a 2-dimensional array.')
  2336. if s[1] != sB[1]:
  2337. raise ValueError('XA and XB must have the same number of columns '
  2338. '(i.e. feature dimension.)')
  2339. mA = s[0]
  2340. mB = sB[0]
  2341. n = s[1]
  2342. if callable(metric):
  2343. mstr = getattr(metric, '__name__', 'Unknown')
  2344. metric_info = _METRIC_ALIAS.get(mstr, None)
  2345. if metric_info is not None:
  2346. XA, XB, typ, kwargs = _validate_cdist_input(
  2347. XA, XB, mA, mB, n, metric_info, **kwargs)
  2348. return _cdist_callable(XA, XB, metric=metric, out=out, **kwargs)
  2349. elif isinstance(metric, str):
  2350. mstr = metric.lower()
  2351. metric_info = _METRIC_ALIAS.get(mstr, None)
  2352. if metric_info is not None:
  2353. cdist_fn = metric_info.cdist_func
  2354. return cdist_fn(XA, XB, out=out, **kwargs)
  2355. elif mstr.startswith("test_"):
  2356. metric_info = _TEST_METRICS.get(mstr, None)
  2357. if metric_info is None:
  2358. raise ValueError(f'Unknown "Test" Distance Metric: {mstr[5:]}')
  2359. XA, XB, typ, kwargs = _validate_cdist_input(
  2360. XA, XB, mA, mB, n, metric_info, **kwargs)
  2361. return _cdist_callable(
  2362. XA, XB, metric=metric_info.dist_func, out=out, **kwargs)
  2363. else:
  2364. raise ValueError('Unknown Distance Metric: %s' % mstr)
  2365. else:
  2366. raise TypeError('2nd argument metric must be a string identifier '
  2367. 'or a function.')