_fitpack_impl.py 28 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824
  1. """
  2. fitpack (dierckx in netlib) --- A Python-C wrapper to FITPACK (by P. Dierckx).
  3. FITPACK is a collection of FORTRAN programs for curve and surface
  4. fitting with splines and tensor product splines.
  5. See
  6. https://web.archive.org/web/20010524124604/http://www.cs.kuleuven.ac.be:80/cwis/research/nalag/research/topics/fitpack.html
  7. or
  8. http://www.netlib.org/dierckx/
  9. Copyright 2002 Pearu Peterson all rights reserved,
  10. Pearu Peterson <pearu@cens.ioc.ee>
  11. Permission to use, modify, and distribute this software is given under the
  12. terms of the SciPy (BSD style) license. See LICENSE.txt that came with
  13. this distribution for specifics.
  14. NO WARRANTY IS EXPRESSED OR IMPLIED. USE AT YOUR OWN RISK.
  15. TODO: Make interfaces to the following fitpack functions:
  16. For univariate splines: cocosp, concon, fourco, insert
  17. For bivariate splines: profil, regrid, parsur, surev
  18. """
  19. __all__ = ['splrep', 'splprep', 'splev', 'splint', 'sproot', 'spalde',
  20. 'bisplrep', 'bisplev', 'insert', 'splder', 'splantider']
  21. import warnings
  22. import numpy as np
  23. from . import _fitpack
  24. from numpy import (atleast_1d, array, ones, zeros, sqrt, ravel, transpose,
  25. empty, iinfo, asarray)
  26. # Try to replace _fitpack interface with
  27. # f2py-generated version
  28. from . import _dfitpack as dfitpack
  29. from scipy._lib._array_api import array_namespace, concat_1d, xp_capabilities
  30. dfitpack_int = dfitpack.types.intvar.dtype
  31. def _int_overflow(x, exception, msg=None):
  32. """Cast the value to an dfitpack_int and raise an OverflowError if the value
  33. cannot fit.
  34. """
  35. if x > iinfo(dfitpack_int).max:
  36. if msg is None:
  37. msg = f'{x!r} cannot fit into an {dfitpack_int!r}'
  38. raise exception(msg)
  39. return dfitpack_int.type(x)
  40. _iermess = {
  41. 0: ["The spline has a residual sum of squares fp such that "
  42. "abs(fp-s)/s<=0.001", None],
  43. -1: ["The spline is an interpolating spline (fp=0)", None],
  44. -2: ["The spline is weighted least-squares polynomial of degree k.\n"
  45. "fp gives the upper bound fp0 for the smoothing factor s", None],
  46. 1: ["The required storage space exceeds the available storage space.\n"
  47. "Probable causes: data (x,y) size is too small or smoothing parameter"
  48. "\ns is too small (fp>s).", ValueError],
  49. 2: ["A theoretically impossible result when finding a smoothing spline\n"
  50. "with fp = s. Probable cause: s too small. (abs(fp-s)/s>0.001)",
  51. ValueError],
  52. 3: ["The maximal number of iterations (20) allowed for finding smoothing\n"
  53. "spline with fp=s has been reached. Probable cause: s too small.\n"
  54. "(abs(fp-s)/s>0.001)", ValueError],
  55. 10: ["Error on input data", ValueError],
  56. 'unknown': ["An error occurred", TypeError]
  57. }
  58. _iermess2 = {
  59. 0: ["The spline has a residual sum of squares fp such that "
  60. "abs(fp-s)/s<=0.001", None],
  61. -1: ["The spline is an interpolating spline (fp=0)", None],
  62. -2: ["The spline is weighted least-squares polynomial of degree kx and ky."
  63. "\nfp gives the upper bound fp0 for the smoothing factor s", None],
  64. -3: ["Warning. The coefficients of the spline have been computed as the\n"
  65. "minimal norm least-squares solution of a rank deficient system.",
  66. None],
  67. 1: ["The required storage space exceeds the available storage space.\n"
  68. "Probable causes: nxest or nyest too small or s is too small. (fp>s)",
  69. ValueError],
  70. 2: ["A theoretically impossible result when finding a smoothing spline\n"
  71. "with fp = s. Probable causes: s too small or badly chosen eps.\n"
  72. "(abs(fp-s)/s>0.001)", ValueError],
  73. 3: ["The maximal number of iterations (20) allowed for finding smoothing\n"
  74. "spline with fp=s has been reached. Probable cause: s too small.\n"
  75. "(abs(fp-s)/s>0.001)", ValueError],
  76. 4: ["No more knots can be added because the number of B-spline\n"
  77. "coefficients already exceeds the number of data points m.\n"
  78. "Probable causes: either s or m too small. (fp>s)", ValueError],
  79. 5: ["No more knots can be added because the additional knot would\n"
  80. "coincide with an old one. Probable cause: s too small or too large\n"
  81. "a weight to an inaccurate data point. (fp>s)", ValueError],
  82. 10: ["Error on input data", ValueError],
  83. 11: ["rwrk2 too small, i.e., there is not enough workspace for computing\n"
  84. "the minimal least-squares solution of a rank deficient system of\n"
  85. "linear equations.", ValueError],
  86. 'unknown': ["An error occurred", TypeError]
  87. }
  88. _parcur_cache = {'t': array([], float), 'wrk': array([], float),
  89. 'iwrk': array([], dfitpack_int), 'u': array([], float),
  90. 'ub': 0, 'ue': 1}
  91. def splprep(x, w=None, u=None, ub=None, ue=None, k=3, task=0, s=None, t=None,
  92. full_output=0, nest=None, per=0, quiet=1):
  93. # see the docstring of `_fitpack_py/splprep`
  94. if task <= 0:
  95. _parcur_cache = {'t': array([], float), 'wrk': array([], float),
  96. 'iwrk': array([], dfitpack_int), 'u': array([], float),
  97. 'ub': 0, 'ue': 1}
  98. x = atleast_1d(x)
  99. idim, m = x.shape
  100. if per:
  101. for i in range(idim):
  102. if x[i][0] != x[i][-1]:
  103. if not quiet:
  104. warnings.warn(
  105. RuntimeWarning(f'Setting x[{i}][{m}]=x[{i}][0]'),
  106. stacklevel=2
  107. )
  108. x[i][-1] = x[i][0]
  109. if not 0 < idim < 11:
  110. raise TypeError('0 < idim < 11 must hold')
  111. if w is None:
  112. w = ones(m, float)
  113. else:
  114. w = atleast_1d(w)
  115. ipar = (u is not None)
  116. if ipar:
  117. _parcur_cache['u'] = u
  118. if ub is None:
  119. _parcur_cache['ub'] = u[0]
  120. else:
  121. _parcur_cache['ub'] = ub
  122. if ue is None:
  123. _parcur_cache['ue'] = u[-1]
  124. else:
  125. _parcur_cache['ue'] = ue
  126. else:
  127. _parcur_cache['u'] = zeros(m, float)
  128. if not (1 <= k <= 5):
  129. raise TypeError(f'1 <= k= {k} <=5 must hold')
  130. if not (-1 <= task <= 1):
  131. raise TypeError('task must be -1, 0 or 1')
  132. if (not len(w) == m) or (ipar == 1 and (not len(u) == m)):
  133. raise TypeError('Mismatch of input dimensions')
  134. if s is None:
  135. s = m - sqrt(2*m)
  136. if t is None and task == -1:
  137. raise TypeError('Knots must be given for task=-1')
  138. if t is not None:
  139. _parcur_cache['t'] = atleast_1d(t)
  140. n = len(_parcur_cache['t'])
  141. if task == -1 and n < 2*k + 2:
  142. raise TypeError('There must be at least 2*k+2 knots for task=-1')
  143. if m <= k:
  144. raise TypeError('m > k must hold')
  145. if nest is None:
  146. nest = m + 2*k
  147. if (task >= 0 and s == 0) or (nest < 0):
  148. if per:
  149. nest = m + 2*k
  150. else:
  151. nest = m + k + 1
  152. nest = max(nest, 2*k + 3)
  153. u = _parcur_cache['u']
  154. ub = _parcur_cache['ub']
  155. ue = _parcur_cache['ue']
  156. t = _parcur_cache['t']
  157. wrk = _parcur_cache['wrk']
  158. iwrk = _parcur_cache['iwrk']
  159. t, c, o = _fitpack._parcur(ravel(transpose(x)), w, u, ub, ue, k,
  160. task, ipar, s, t, nest, wrk, iwrk, per)
  161. _parcur_cache['u'] = o['u']
  162. _parcur_cache['ub'] = o['ub']
  163. _parcur_cache['ue'] = o['ue']
  164. _parcur_cache['t'] = t
  165. _parcur_cache['wrk'] = o['wrk']
  166. _parcur_cache['iwrk'] = o['iwrk']
  167. ier = o['ier']
  168. fp = o['fp']
  169. n = len(t)
  170. u = o['u']
  171. c = c.reshape((idim, n - k - 1))
  172. tcku = [t, list(c), k], u
  173. if ier <= 0 and not quiet:
  174. warnings.warn(
  175. RuntimeWarning(
  176. _iermess[ier][0] + f"\tk={k} n={len(t)} m={m} fp={fp} s={s}"
  177. ),
  178. stacklevel=2
  179. )
  180. if ier > 0 and not full_output:
  181. if ier in [1, 2, 3]:
  182. warnings.warn(RuntimeWarning(_iermess[ier][0]), stacklevel=2)
  183. else:
  184. try:
  185. raise _iermess[ier][1](_iermess[ier][0])
  186. except KeyError as e:
  187. raise _iermess['unknown'][1](_iermess['unknown'][0]) from e
  188. if full_output:
  189. try:
  190. return tcku, fp, ier, _iermess[ier][0]
  191. except KeyError:
  192. return tcku, fp, ier, _iermess['unknown'][0]
  193. else:
  194. return tcku
  195. _curfit_cache = {'t': array([], float), 'wrk': array([], float),
  196. 'iwrk': array([], dfitpack_int)}
  197. def splrep(x, y, w=None, xb=None, xe=None, k=3, task=0, s=None, t=None,
  198. full_output=0, per=0, quiet=1):
  199. # see the docstring of `_fitpack_py/splrep`
  200. if task <= 0:
  201. _curfit_cache = {}
  202. x, y = map(atleast_1d, [x, y])
  203. m = len(x)
  204. if w is None:
  205. w = ones(m, float)
  206. if s is None:
  207. s = 0.0
  208. else:
  209. w = atleast_1d(w)
  210. if s is None:
  211. s = m - sqrt(2*m)
  212. if not len(w) == m:
  213. raise TypeError(f'len(w)={len(w)} is not equal to m={m}')
  214. if (m != len(y)) or (m != len(w)):
  215. raise TypeError('Lengths of the first three arguments (x,y,w) must '
  216. 'be equal')
  217. if not (1 <= k <= 5):
  218. raise TypeError(
  219. f'Given degree of the spline (k={k}) is not supported. (1<=k<=5)'
  220. )
  221. if m <= k:
  222. raise TypeError('m > k must hold')
  223. if xb is None:
  224. xb = x[0]
  225. if xe is None:
  226. xe = x[-1]
  227. if not (-1 <= task <= 1):
  228. raise TypeError('task must be -1, 0 or 1')
  229. if t is not None:
  230. task = -1
  231. if task == -1:
  232. if t is None:
  233. raise TypeError('Knots must be given for task=-1')
  234. numknots = len(t)
  235. _curfit_cache['t'] = empty((numknots + 2*k + 2,), float)
  236. _curfit_cache['t'][k+1:-k-1] = t
  237. nest = len(_curfit_cache['t'])
  238. elif task == 0:
  239. if per:
  240. nest = max(m + 2*k, 2*k + 3)
  241. else:
  242. nest = max(m + k + 1, 2*k + 3)
  243. t = empty((nest,), float)
  244. _curfit_cache['t'] = t
  245. if task <= 0:
  246. if per:
  247. _curfit_cache['wrk'] = empty((m*(k + 1) + nest*(8 + 5*k),), float)
  248. else:
  249. _curfit_cache['wrk'] = empty((m*(k + 1) + nest*(7 + 3*k),), float)
  250. _curfit_cache['iwrk'] = empty((nest,), dfitpack_int)
  251. try:
  252. t = _curfit_cache['t']
  253. wrk = _curfit_cache['wrk']
  254. iwrk = _curfit_cache['iwrk']
  255. except KeyError as e:
  256. raise TypeError("must call with task=1 only after"
  257. " call with task=0,-1") from e
  258. if not per:
  259. n, c, fp, ier = dfitpack.curfit(task, x, y, w, t, wrk, iwrk,
  260. xb, xe, k, s)
  261. else:
  262. n, c, fp, ier = dfitpack.percur(task, x, y, w, t, wrk, iwrk, k, s)
  263. tck = (t[:n], c[:n], k)
  264. if ier <= 0 and not quiet:
  265. _mess = (_iermess[ier][0] + f"\tk={k} n={len(t)} m={m} fp={fp} s={s}")
  266. warnings.warn(RuntimeWarning(_mess), stacklevel=2)
  267. if ier > 0 and not full_output:
  268. if ier in [1, 2, 3]:
  269. warnings.warn(RuntimeWarning(_iermess[ier][0]), stacklevel=2)
  270. else:
  271. try:
  272. raise _iermess[ier][1](_iermess[ier][0])
  273. except KeyError as e:
  274. raise _iermess['unknown'][1](_iermess['unknown'][0]) from e
  275. if full_output:
  276. try:
  277. return tck, fp, ier, _iermess[ier][0]
  278. except KeyError:
  279. return tck, fp, ier, _iermess['unknown'][0]
  280. else:
  281. return tck
  282. def splev(x, tck, der=0, ext=0):
  283. # see the docstring of `_fitpack_py/splev`
  284. t, c, k = tck
  285. try:
  286. c[0][0]
  287. parametric = True
  288. except Exception:
  289. parametric = False
  290. if parametric:
  291. return list(map(lambda c, x=x, t=t, k=k, der=der:
  292. splev(x, [t, c, k], der, ext), c))
  293. else:
  294. if not (0 <= der <= k):
  295. raise ValueError(f"0<=der={der}<=k={k} must hold")
  296. if ext not in (0, 1, 2, 3):
  297. raise ValueError(f"ext = {ext} not in (0, 1, 2, 3) ")
  298. x = asarray(x)
  299. shape = x.shape
  300. x = atleast_1d(x).ravel()
  301. if der == 0:
  302. y, ier = dfitpack.splev(t, c, k, x, ext)
  303. else:
  304. y, ier = dfitpack.splder(t, c, k, x, der, ext)
  305. if ier == 10:
  306. raise ValueError("Invalid input data")
  307. if ier == 1:
  308. raise ValueError("Found x value not in the domain")
  309. if ier:
  310. raise TypeError("An error occurred")
  311. return y.reshape(shape)
  312. def splint(a, b, tck, full_output=0):
  313. # see the docstring of `_fitpack_py/splint`
  314. t, c, k = tck
  315. try:
  316. c[0][0]
  317. parametric = True
  318. except Exception:
  319. parametric = False
  320. if parametric:
  321. return list(map(lambda c, a=a, b=b, t=t, k=k:
  322. splint(a, b, [t, c, k]), c))
  323. else:
  324. aint, wrk = dfitpack.splint(t, c, k, a, b)
  325. if full_output:
  326. return aint, wrk
  327. else:
  328. return aint
  329. def sproot(tck, mest=10):
  330. # see the docstring of `_fitpack_py/sproot`
  331. t, c, k = tck
  332. if k != 3:
  333. raise ValueError("sproot works only for cubic (k=3) splines")
  334. try:
  335. c[0][0]
  336. parametric = True
  337. except Exception:
  338. parametric = False
  339. if parametric:
  340. return list(map(lambda c, t=t, k=k, mest=mest:
  341. sproot([t, c, k], mest), c))
  342. else:
  343. if len(t) < 8:
  344. raise TypeError(f"The number of knots {len(t)}>=8")
  345. z, m, ier = dfitpack.sproot(t, c, mest)
  346. if ier == 10:
  347. raise TypeError("Invalid input data. "
  348. "t1<=..<=t4<t5<..<tn-3<=..<=tn must hold.")
  349. if ier == 0:
  350. return z[:m]
  351. if ier == 1:
  352. warnings.warn(RuntimeWarning("The number of zeros exceeds mest"),
  353. stacklevel=2)
  354. return z[:m]
  355. raise TypeError("Unknown error")
  356. def spalde(x, tck):
  357. # see the docstring of `_fitpack_py/spalde`
  358. t, c, k = tck
  359. try:
  360. c[0][0]
  361. parametric = True
  362. except Exception:
  363. parametric = False
  364. if parametric:
  365. return list(map(lambda c, x=x, t=t, k=k:
  366. spalde(x, [t, c, k]), c))
  367. else:
  368. x = atleast_1d(x)
  369. if len(x) > 1:
  370. return list(map(lambda x, tck=tck: spalde(x, tck), x))
  371. d, ier = dfitpack.spalde(t, c, k+1, x[0])
  372. if ier == 0:
  373. return d
  374. if ier == 10:
  375. raise TypeError("Invalid input data. t(k)<=x<=t(n-k+1) must hold.")
  376. raise TypeError("Unknown error")
  377. # def _curfit(x,y,w=None,xb=None,xe=None,k=3,task=0,s=None,t=None,
  378. # full_output=0,nest=None,per=0,quiet=1):
  379. _surfit_cache = {'tx': array([], float), 'ty': array([], float),
  380. 'wrk': array([], float), 'iwrk': array([], dfitpack_int)}
  381. @xp_capabilities(out_of_scope=True)
  382. def bisplrep(x, y, z, w=None, xb=None, xe=None, yb=None, ye=None,
  383. kx=3, ky=3, task=0, s=None, eps=1e-16, tx=None, ty=None,
  384. full_output=0, nxest=None, nyest=None, quiet=1):
  385. """
  386. Find a bivariate B-spline representation of a surface.
  387. Given a set of data points (x[i], y[i], z[i]) representing a surface
  388. z=f(x,y), compute a B-spline representation of the surface. Based on
  389. the routine SURFIT from FITPACK.
  390. Parameters
  391. ----------
  392. x, y, z : ndarray
  393. Rank-1 arrays of data points.
  394. w : ndarray, optional
  395. Rank-1 array of weights. By default ``w=np.ones(len(x))``.
  396. xb, xe : float, optional
  397. End points of approximation interval in `x`.
  398. By default ``xb = x.min(), xe=x.max()``.
  399. yb, ye : float, optional
  400. End points of approximation interval in `y`.
  401. By default ``yb=y.min(), ye = y.max()``.
  402. kx, ky : int, optional
  403. The degrees of the spline (1 <= kx, ky <= 5).
  404. Third order (kx=ky=3) is recommended.
  405. task : int, optional
  406. If task=0, find knots in x and y and coefficients for a given
  407. smoothing factor, s.
  408. If task=1, find knots and coefficients for another value of the
  409. smoothing factor, s. bisplrep must have been previously called
  410. with task=0 or task=1.
  411. If task=-1, find coefficients for a given set of knots tx, ty.
  412. s : float, optional
  413. A non-negative smoothing factor. If weights correspond
  414. to the inverse of the standard-deviation of the errors in z,
  415. then a good s-value should be found in the range
  416. ``(m-sqrt(2*m),m+sqrt(2*m))`` where m=len(x).
  417. eps : float, optional
  418. A threshold for determining the effective rank of an
  419. over-determined linear system of equations (0 < eps < 1).
  420. `eps` is not likely to need changing.
  421. tx, ty : ndarray, optional
  422. Rank-1 arrays of the knots of the spline for task=-1
  423. full_output : int, optional
  424. Non-zero to return optional outputs.
  425. nxest, nyest : int, optional
  426. Over-estimates of the total number of knots. If None then
  427. ``nxest = max(kx+sqrt(m/2),2*kx+3)``,
  428. ``nyest = max(ky+sqrt(m/2),2*ky+3)``.
  429. quiet : int, optional
  430. Non-zero to suppress printing of messages.
  431. Returns
  432. -------
  433. tck : array_like
  434. A list [tx, ty, c, kx, ky] containing the knots (tx, ty) and
  435. coefficients (c) of the bivariate B-spline representation of the
  436. surface along with the degree of the spline.
  437. fp : ndarray
  438. The weighted sum of squared residuals of the spline approximation.
  439. ier : int
  440. An integer flag about splrep success. Success is indicated if
  441. ier<=0. If ier in [1,2,3] an error occurred but was not raised.
  442. Otherwise an error is raised.
  443. msg : str
  444. A message corresponding to the integer flag, ier.
  445. See Also
  446. --------
  447. splprep, splrep, splint, sproot, splev
  448. UnivariateSpline, BivariateSpline
  449. Notes
  450. -----
  451. See `bisplev` to evaluate the value of the B-spline given its tck
  452. representation.
  453. If the input data is such that input dimensions have incommensurate
  454. units and differ by many orders of magnitude, the interpolant may have
  455. numerical artifacts. Consider rescaling the data before interpolation.
  456. References
  457. ----------
  458. .. [1] Dierckx P.:An algorithm for surface fitting with spline functions
  459. Ima J. Numer. Anal. 1 (1981) 267-283.
  460. .. [2] Dierckx P.:An algorithm for surface fitting with spline functions
  461. report tw50, Dept. Computer Science,K.U.Leuven, 1980.
  462. .. [3] Dierckx P.:Curve and surface fitting with splines, Monographs on
  463. Numerical Analysis, Oxford University Press, 1993.
  464. Examples
  465. --------
  466. Examples are given :ref:`in the tutorial <tutorial-interpolate_2d_spline>`.
  467. """
  468. x, y, z = map(ravel, [x, y, z]) # ensure 1-d arrays.
  469. m = len(x)
  470. if not (m == len(y) == len(z)):
  471. raise TypeError('len(x)==len(y)==len(z) must hold.')
  472. if w is None:
  473. w = ones(m, float)
  474. else:
  475. w = atleast_1d(w)
  476. if not len(w) == m:
  477. raise TypeError(f'len(w)={len(w)} is not equal to m={m}')
  478. if xb is None:
  479. xb = x.min()
  480. if xe is None:
  481. xe = x.max()
  482. if yb is None:
  483. yb = y.min()
  484. if ye is None:
  485. ye = y.max()
  486. if not (-1 <= task <= 1):
  487. raise TypeError('task must be -1, 0 or 1')
  488. if s is None:
  489. s = m - sqrt(2*m)
  490. if tx is None and task == -1:
  491. raise TypeError('Knots_x must be given for task=-1')
  492. if tx is not None:
  493. _surfit_cache['tx'] = atleast_1d(tx)
  494. nx = len(_surfit_cache['tx'])
  495. if ty is None and task == -1:
  496. raise TypeError('Knots_y must be given for task=-1')
  497. if ty is not None:
  498. _surfit_cache['ty'] = atleast_1d(ty)
  499. ny = len(_surfit_cache['ty'])
  500. if task == -1 and nx < 2*kx+2:
  501. raise TypeError('There must be at least 2*kx+2 knots_x for task=-1')
  502. if task == -1 and ny < 2*ky+2:
  503. raise TypeError('There must be at least 2*ky+2 knots_x for task=-1')
  504. if not ((1 <= kx <= 5) and (1 <= ky <= 5)):
  505. raise TypeError(
  506. f'Given degree of the spline (kx,ky={kx},{ky}) is not supported. (1<=k<=5)'
  507. )
  508. if m < (kx + 1)*(ky + 1):
  509. raise TypeError('m >= (kx+1)(ky+1) must hold')
  510. if nxest is None:
  511. nxest = int(kx + sqrt(m/2))
  512. if nyest is None:
  513. nyest = int(ky + sqrt(m/2))
  514. nxest, nyest = max(nxest, 2*kx + 3), max(nyest, 2*ky + 3)
  515. if task >= 0 and s == 0:
  516. nxest = int(kx + sqrt(3*m))
  517. nyest = int(ky + sqrt(3*m))
  518. if task == -1:
  519. _surfit_cache['tx'] = atleast_1d(tx)
  520. _surfit_cache['ty'] = atleast_1d(ty)
  521. tx, ty = _surfit_cache['tx'], _surfit_cache['ty']
  522. wrk = _surfit_cache['wrk']
  523. u = nxest - kx - 1
  524. v = nyest - ky - 1
  525. km = max(kx, ky) + 1
  526. ne = max(nxest, nyest)
  527. bx, by = kx*v + ky + 1, ky*u + kx + 1
  528. b1, b2 = bx, bx + v - ky
  529. if bx > by:
  530. b1, b2 = by, by + u - kx
  531. msg = "Too many data points to interpolate"
  532. lwrk1 = _int_overflow(u*v*(2 + b1 + b2) +
  533. 2*(u + v + km*(m + ne) + ne - kx - ky) + b2 + 1,
  534. OverflowError,
  535. msg=msg)
  536. lwrk2 = _int_overflow(u*v*(b2 + 1) + b2, OverflowError, msg=msg)
  537. tx, ty, c, o = _fitpack._surfit(x, y, z, w, xb, xe, yb, ye, kx, ky,
  538. task, s, eps, tx, ty, nxest, nyest,
  539. wrk, lwrk1, lwrk2)
  540. _curfit_cache['tx'] = tx
  541. _curfit_cache['ty'] = ty
  542. _curfit_cache['wrk'] = o['wrk']
  543. ier, fp = o['ier'], o['fp']
  544. tck = [tx, ty, c, kx, ky]
  545. ierm = min(11, max(-3, ier))
  546. if ierm <= 0 and not quiet:
  547. _mess = (
  548. _iermess2[ierm][0] +
  549. f"\tkx,ky={kx},{ky} nx,ny={len(tx)},{len(ty)} m={m} fp={fp} s={s}"
  550. )
  551. warnings.warn(RuntimeWarning(_mess), stacklevel=2)
  552. if ierm > 0 and not full_output:
  553. if ier in [1, 2, 3, 4, 5]:
  554. _mess = (
  555. f"\n\tkx,ky={kx},{ky} nx,ny={len(tx)},{len(ty)} m={m} fp={fp} s={s}"
  556. )
  557. warnings.warn(RuntimeWarning(_iermess2[ierm][0] + _mess), stacklevel=2)
  558. else:
  559. try:
  560. raise _iermess2[ierm][1](_iermess2[ierm][0])
  561. except KeyError as e:
  562. raise _iermess2['unknown'][1](_iermess2['unknown'][0]) from e
  563. if full_output:
  564. try:
  565. return tck, fp, ier, _iermess2[ierm][0]
  566. except KeyError:
  567. return tck, fp, ier, _iermess2['unknown'][0]
  568. else:
  569. return tck
  570. @xp_capabilities(out_of_scope=True)
  571. def bisplev(x, y, tck, dx=0, dy=0):
  572. """
  573. Evaluate a bivariate B-spline and its derivatives.
  574. Return a rank-2 array of spline function values (or spline derivative
  575. values) at points given by the cross-product of the rank-1 arrays `x` and
  576. `y`. In special cases, return an array or just a float if either `x` or
  577. `y` or both are floats. Based on BISPEV and PARDER from FITPACK.
  578. Parameters
  579. ----------
  580. x, y : ndarray
  581. Rank-1 arrays specifying the domain over which to evaluate the
  582. spline or its derivative.
  583. tck : tuple
  584. A sequence of length 5 returned by `bisplrep` containing the knot
  585. locations, the coefficients, and the degree of the spline:
  586. [tx, ty, c, kx, ky].
  587. dx, dy : int, optional
  588. The orders of the partial derivatives in `x` and `y` respectively.
  589. Returns
  590. -------
  591. vals : ndarray
  592. The B-spline or its derivative evaluated over the set formed by
  593. the cross-product of `x` and `y`.
  594. See Also
  595. --------
  596. splprep, splrep, splint, sproot, splev
  597. UnivariateSpline, BivariateSpline
  598. Notes
  599. -----
  600. See `bisplrep` to generate the `tck` representation.
  601. References
  602. ----------
  603. .. [1] Dierckx P. : An algorithm for surface fitting
  604. with spline functions
  605. Ima J. Numer. Anal. 1 (1981) 267-283.
  606. .. [2] Dierckx P. : An algorithm for surface fitting
  607. with spline functions
  608. report tw50, Dept. Computer Science,K.U.Leuven, 1980.
  609. .. [3] Dierckx P. : Curve and surface fitting with splines,
  610. Monographs on Numerical Analysis, Oxford University Press, 1993.
  611. Examples
  612. --------
  613. Examples are given :ref:`in the tutorial <tutorial-interpolate_2d_spline>`.
  614. """
  615. tx, ty, c, kx, ky = tck
  616. if not (0 <= dx < kx):
  617. raise ValueError(f"0 <= dx = {dx} < kx = {kx} must hold")
  618. if not (0 <= dy < ky):
  619. raise ValueError(f"0 <= dy = {dy} < ky = {ky} must hold")
  620. x, y = map(atleast_1d, [x, y])
  621. if (len(x.shape) != 1) or (len(y.shape) != 1):
  622. raise ValueError("First two entries should be rank-1 arrays.")
  623. msg = "Too many data points to interpolate."
  624. _int_overflow(x.size * y.size, MemoryError, msg=msg)
  625. if dx != 0 or dy != 0:
  626. _int_overflow((tx.size - kx - 1)*(ty.size - ky - 1),
  627. MemoryError, msg=msg)
  628. z, ier = dfitpack.parder(tx, ty, c, kx, ky, dx, dy, x, y)
  629. else:
  630. z, ier = dfitpack.bispev(tx, ty, c, kx, ky, x, y)
  631. if ier == 10:
  632. raise ValueError("Invalid input data")
  633. if ier:
  634. raise TypeError("An error occurred")
  635. z = z.reshape((len(x), len(y)))
  636. if len(z) > 1:
  637. return z
  638. if len(z[0]) > 1:
  639. return z[0]
  640. return z[0][0]
  641. def dblint(xa, xb, ya, yb, tck):
  642. """Evaluate the integral of a spline over area [xa,xb] x [ya,yb].
  643. Parameters
  644. ----------
  645. xa, xb : float
  646. The end-points of the x integration interval.
  647. ya, yb : float
  648. The end-points of the y integration interval.
  649. tck : list [tx, ty, c, kx, ky]
  650. A sequence of length 5 returned by bisplrep containing the knot
  651. locations tx, ty, the coefficients c, and the degrees kx, ky
  652. of the spline.
  653. Returns
  654. -------
  655. integ : float
  656. The value of the resulting integral.
  657. """
  658. tx, ty, c, kx, ky = tck
  659. return dfitpack.dblint(tx, ty, c, kx, ky, xa, xb, ya, yb)
  660. def insert(x, tck, m=1, per=0):
  661. # see the docstring of `_fitpack_py/insert`
  662. t, c, k = tck
  663. try:
  664. c[0][0]
  665. parametric = True
  666. except Exception:
  667. parametric = False
  668. if parametric:
  669. cc = []
  670. for c_vals in c:
  671. tt, cc_val, kk = insert(x, [t, c_vals, k], m)
  672. cc.append(cc_val)
  673. return (tt, cc, kk)
  674. else:
  675. tt, cc, ier = _fitpack._insert(per, t, c, k, x, m)
  676. if ier == 10:
  677. raise ValueError("Invalid input data")
  678. if ier:
  679. raise TypeError("An error occurred")
  680. return (tt, cc, k)
  681. def splder(tck, n=1, xp=None):
  682. # see the docstring of `_fitpack_py/splder`
  683. if n < 0:
  684. return splantider(tck, -n)
  685. t, c, k = tck
  686. if xp is None:
  687. xp = array_namespace(t, c)
  688. if n > k:
  689. raise ValueError(f"Order of derivative (n = {n!r}) must be <= "
  690. f"order of spline (k = {tck[2]!r})")
  691. # Extra axes for the trailing dims of the `c` array:
  692. sh = (slice(None),) + ((None,)*len(c.shape[1:]))
  693. with np.errstate(invalid='raise', divide='raise'):
  694. try:
  695. for j in range(n):
  696. # See e.g. Schumaker, Spline Functions: Basic Theory, Chapter 5
  697. # Compute the denominator in the differentiation formula.
  698. # (and append trailing dims, if necessary)
  699. dt = t[k+1:-1] - t[1:-k-1]
  700. dt = dt[sh]
  701. # Compute the new coefficients
  702. c = (c[1:-1-k, ...] - c[:-2-k, ...]) * k / dt
  703. # Pad coefficient array to same size as knots (FITPACK
  704. # convention)
  705. c = concat_1d(xp, c, xp.zeros((k,) + c.shape[1:]))
  706. # Adjust knots
  707. t = t[1:-1]
  708. k -= 1
  709. except FloatingPointError as e:
  710. raise ValueError("The spline has internal repeated knots "
  711. f"and is not differentiable {n} times") from e
  712. return t, c, k
  713. def splantider(tck, n=1, *, xp=None):
  714. # see the docstring of `_fitpack_py/splantider`
  715. if n < 0:
  716. return splder(tck, -n)
  717. t, c, k = tck
  718. if xp is None:
  719. xp = array_namespace(t, c)
  720. # Extra axes for the trailing dims of the `c` array:
  721. sh = (slice(None),) + (None,)*len(c.shape[1:])
  722. for j in range(n):
  723. # This is the inverse set of operations to splder.
  724. # Compute the multiplier in the antiderivative formula.
  725. dt = t[k+1:] - t[:-k-1]
  726. dt = dt[sh]
  727. # Compute the new coefficients
  728. c = xp.cumulative_sum(c[:-k-1, ...] * dt, axis=0) / (k + 1)
  729. c = concat_1d(
  730. xp,
  731. xp.zeros((1,) + c.shape[1:]),
  732. c,
  733. xp.stack([c[-1, ...]] * (k+2)),
  734. )
  735. # New knots
  736. t = concat_1d(xp, t[0], t, t[-1])
  737. k += 1
  738. return t, c, k