fitpack2.py 61 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716
  1. """
  2. fitpack --- curve and surface fitting with splines
  3. fitpack is based on a collection of Fortran routines DIERCKX
  4. by P. Dierckx (see http://www.netlib.org/dierckx/) transformed
  5. to double routines by Pearu Peterson.
  6. """
  7. # Created by Pearu Peterson, June,August 2003
  8. from __future__ import division, print_function, absolute_import
  9. __all__ = [
  10. 'UnivariateSpline',
  11. 'InterpolatedUnivariateSpline',
  12. 'LSQUnivariateSpline',
  13. 'BivariateSpline',
  14. 'LSQBivariateSpline',
  15. 'SmoothBivariateSpline',
  16. 'LSQSphereBivariateSpline',
  17. 'SmoothSphereBivariateSpline',
  18. 'RectBivariateSpline',
  19. 'RectSphereBivariateSpline']
  20. import warnings
  21. from numpy import zeros, concatenate, alltrue, ravel, all, diff, array, ones
  22. import numpy as np
  23. from . import fitpack
  24. from . import dfitpack
  25. # ############### Univariate spline ####################
  26. _curfit_messages = {1: """
  27. The required storage space exceeds the available storage space, as
  28. specified by the parameter nest: nest too small. If nest is already
  29. large (say nest > m/2), it may also indicate that s is too small.
  30. The approximation returned is the weighted least-squares spline
  31. according to the knots t[0],t[1],...,t[n-1]. (n=nest) the parameter fp
  32. gives the corresponding weighted sum of squared residuals (fp>s).
  33. """,
  34. 2: """
  35. A theoretically impossible result was found during the iteration
  36. process for finding a smoothing spline with fp = s: s too small.
  37. There is an approximation returned but the corresponding weighted sum
  38. of squared residuals does not satisfy the condition abs(fp-s)/s < tol.""",
  39. 3: """
  40. The maximal number of iterations maxit (set to 20 by the program)
  41. allowed for finding a smoothing spline with fp=s has been reached: s
  42. too small.
  43. There is an approximation returned but the corresponding weighted sum
  44. of squared residuals does not satisfy the condition abs(fp-s)/s < tol.""",
  45. 10: """
  46. Error on entry, no approximation returned. The following conditions
  47. must hold:
  48. xb<=x[0]<x[1]<...<x[m-1]<=xe, w[i]>0, i=0..m-1
  49. if iopt=-1:
  50. xb<t[k+1]<t[k+2]<...<t[n-k-2]<xe"""
  51. }
  52. # UnivariateSpline, ext parameter can be an int or a string
  53. _extrap_modes = {0: 0, 'extrapolate': 0,
  54. 1: 1, 'zeros': 1,
  55. 2: 2, 'raise': 2,
  56. 3: 3, 'const': 3}
  57. class UnivariateSpline(object):
  58. """
  59. One-dimensional smoothing spline fit to a given set of data points.
  60. Fits a spline y = spl(x) of degree `k` to the provided `x`, `y` data. `s`
  61. specifies the number of knots by specifying a smoothing condition.
  62. Parameters
  63. ----------
  64. x : (N,) array_like
  65. 1-D array of independent input data. Must be increasing.
  66. y : (N,) array_like
  67. 1-D array of dependent input data, of the same length as `x`.
  68. w : (N,) array_like, optional
  69. Weights for spline fitting. Must be positive. If None (default),
  70. weights are all equal.
  71. bbox : (2,) array_like, optional
  72. 2-sequence specifying the boundary of the approximation interval. If
  73. None (default), ``bbox=[x[0], x[-1]]``.
  74. k : int, optional
  75. Degree of the smoothing spline. Must be <= 5.
  76. Default is k=3, a cubic spline.
  77. s : float or None, optional
  78. Positive smoothing factor used to choose the number of knots. Number
  79. of knots will be increased until the smoothing condition is satisfied::
  80. sum((w[i] * (y[i]-spl(x[i])))**2, axis=0) <= s
  81. If None (default), ``s = len(w)`` which should be a good value if
  82. ``1/w[i]`` is an estimate of the standard deviation of ``y[i]``.
  83. If 0, spline will interpolate through all data points.
  84. ext : int or str, optional
  85. Controls the extrapolation mode for elements
  86. not in the interval defined by the knot sequence.
  87. * if ext=0 or 'extrapolate', return the extrapolated value.
  88. * if ext=1 or 'zeros', return 0
  89. * if ext=2 or 'raise', raise a ValueError
  90. * if ext=3 of 'const', return the boundary value.
  91. The default value is 0.
  92. check_finite : bool, optional
  93. Whether to check that the input arrays contain only finite numbers.
  94. Disabling may give a performance gain, but may result in problems
  95. (crashes, non-termination or non-sensical results) if the inputs
  96. do contain infinities or NaNs.
  97. Default is False.
  98. See Also
  99. --------
  100. InterpolatedUnivariateSpline : Subclass with smoothing forced to 0
  101. LSQUnivariateSpline : Subclass in which knots are user-selected instead of
  102. being set by smoothing condition
  103. splrep : An older, non object-oriented wrapping of FITPACK
  104. splev, sproot, splint, spalde
  105. BivariateSpline : A similar class for two-dimensional spline interpolation
  106. Notes
  107. -----
  108. The number of data points must be larger than the spline degree `k`.
  109. **NaN handling**: If the input arrays contain ``nan`` values, the result
  110. is not useful, since the underlying spline fitting routines cannot deal
  111. with ``nan`` . A workaround is to use zero weights for not-a-number
  112. data points:
  113. >>> from scipy.interpolate import UnivariateSpline
  114. >>> x, y = np.array([1, 2, 3, 4]), np.array([1, np.nan, 3, 4])
  115. >>> w = np.isnan(y)
  116. >>> y[w] = 0.
  117. >>> spl = UnivariateSpline(x, y, w=~w)
  118. Notice the need to replace a ``nan`` by a numerical value (precise value
  119. does not matter as long as the corresponding weight is zero.)
  120. Examples
  121. --------
  122. >>> import matplotlib.pyplot as plt
  123. >>> from scipy.interpolate import UnivariateSpline
  124. >>> x = np.linspace(-3, 3, 50)
  125. >>> y = np.exp(-x**2) + 0.1 * np.random.randn(50)
  126. >>> plt.plot(x, y, 'ro', ms=5)
  127. Use the default value for the smoothing parameter:
  128. >>> spl = UnivariateSpline(x, y)
  129. >>> xs = np.linspace(-3, 3, 1000)
  130. >>> plt.plot(xs, spl(xs), 'g', lw=3)
  131. Manually change the amount of smoothing:
  132. >>> spl.set_smoothing_factor(0.5)
  133. >>> plt.plot(xs, spl(xs), 'b', lw=3)
  134. >>> plt.show()
  135. """
  136. def __init__(self, x, y, w=None, bbox=[None]*2, k=3, s=None,
  137. ext=0, check_finite=False):
  138. if check_finite:
  139. w_finite = np.isfinite(w).all() if w is not None else True
  140. if (not np.isfinite(x).all() or not np.isfinite(y).all() or
  141. not w_finite):
  142. raise ValueError("x and y array must not contain "
  143. "NaNs or infs.")
  144. if not all(diff(x) > 0.0):
  145. raise ValueError('x must be strictly increasing')
  146. # _data == x,y,w,xb,xe,k,s,n,t,c,fp,fpint,nrdata,ier
  147. try:
  148. self.ext = _extrap_modes[ext]
  149. except KeyError:
  150. raise ValueError("Unknown extrapolation mode %s." % ext)
  151. data = dfitpack.fpcurf0(x, y, k, w=w, xb=bbox[0],
  152. xe=bbox[1], s=s)
  153. if data[-1] == 1:
  154. # nest too small, setting to maximum bound
  155. data = self._reset_nest(data)
  156. self._data = data
  157. self._reset_class()
  158. @classmethod
  159. def _from_tck(cls, tck, ext=0):
  160. """Construct a spline object from given tck"""
  161. self = cls.__new__(cls)
  162. t, c, k = tck
  163. self._eval_args = tck
  164. # _data == x,y,w,xb,xe,k,s,n,t,c,fp,fpint,nrdata,ier
  165. self._data = (None, None, None, None, None, k, None, len(t), t,
  166. c, None, None, None, None)
  167. self.ext = ext
  168. return self
  169. def _reset_class(self):
  170. data = self._data
  171. n, t, c, k, ier = data[7], data[8], data[9], data[5], data[-1]
  172. self._eval_args = t[:n], c[:n], k
  173. if ier == 0:
  174. # the spline returned has a residual sum of squares fp
  175. # such that abs(fp-s)/s <= tol with tol a relative
  176. # tolerance set to 0.001 by the program
  177. pass
  178. elif ier == -1:
  179. # the spline returned is an interpolating spline
  180. self._set_class(InterpolatedUnivariateSpline)
  181. elif ier == -2:
  182. # the spline returned is the weighted least-squares
  183. # polynomial of degree k. In this extreme case fp gives
  184. # the upper bound fp0 for the smoothing factor s.
  185. self._set_class(LSQUnivariateSpline)
  186. else:
  187. # error
  188. if ier == 1:
  189. self._set_class(LSQUnivariateSpline)
  190. message = _curfit_messages.get(ier, 'ier=%s' % (ier))
  191. warnings.warn(message)
  192. def _set_class(self, cls):
  193. self._spline_class = cls
  194. if self.__class__ in (UnivariateSpline, InterpolatedUnivariateSpline,
  195. LSQUnivariateSpline):
  196. self.__class__ = cls
  197. else:
  198. # It's an unknown subclass -- don't change class. cf. #731
  199. pass
  200. def _reset_nest(self, data, nest=None):
  201. n = data[10]
  202. if nest is None:
  203. k, m = data[5], len(data[0])
  204. nest = m+k+1 # this is the maximum bound for nest
  205. else:
  206. if not n <= nest:
  207. raise ValueError("`nest` can only be increased")
  208. t, c, fpint, nrdata = [np.resize(data[j], nest) for j in
  209. [8, 9, 11, 12]]
  210. args = data[:8] + (t, c, n, fpint, nrdata, data[13])
  211. data = dfitpack.fpcurf1(*args)
  212. return data
  213. def set_smoothing_factor(self, s):
  214. """ Continue spline computation with the given smoothing
  215. factor s and with the knots found at the last call.
  216. This routine modifies the spline in place.
  217. """
  218. data = self._data
  219. if data[6] == -1:
  220. warnings.warn('smoothing factor unchanged for'
  221. 'LSQ spline with fixed knots')
  222. return
  223. args = data[:6] + (s,) + data[7:]
  224. data = dfitpack.fpcurf1(*args)
  225. if data[-1] == 1:
  226. # nest too small, setting to maximum bound
  227. data = self._reset_nest(data)
  228. self._data = data
  229. self._reset_class()
  230. def __call__(self, x, nu=0, ext=None):
  231. """
  232. Evaluate spline (or its nu-th derivative) at positions x.
  233. Parameters
  234. ----------
  235. x : array_like
  236. A 1-D array of points at which to return the value of the smoothed
  237. spline or its derivatives. Note: x can be unordered but the
  238. evaluation is more efficient if x is (partially) ordered.
  239. nu : int
  240. The order of derivative of the spline to compute.
  241. ext : int
  242. Controls the value returned for elements of ``x`` not in the
  243. interval defined by the knot sequence.
  244. * if ext=0 or 'extrapolate', return the extrapolated value.
  245. * if ext=1 or 'zeros', return 0
  246. * if ext=2 or 'raise', raise a ValueError
  247. * if ext=3 or 'const', return the boundary value.
  248. The default value is 0, passed from the initialization of
  249. UnivariateSpline.
  250. """
  251. x = np.asarray(x)
  252. # empty input yields empty output
  253. if x.size == 0:
  254. return array([])
  255. # if nu is None:
  256. # return dfitpack.splev(*(self._eval_args+(x,)))
  257. # return dfitpack.splder(nu=nu,*(self._eval_args+(x,)))
  258. if ext is None:
  259. ext = self.ext
  260. else:
  261. try:
  262. ext = _extrap_modes[ext]
  263. except KeyError:
  264. raise ValueError("Unknown extrapolation mode %s." % ext)
  265. return fitpack.splev(x, self._eval_args, der=nu, ext=ext)
  266. def get_knots(self):
  267. """ Return positions of interior knots of the spline.
  268. Internally, the knot vector contains ``2*k`` additional boundary knots.
  269. """
  270. data = self._data
  271. k, n = data[5], data[7]
  272. return data[8][k:n-k]
  273. def get_coeffs(self):
  274. """Return spline coefficients."""
  275. data = self._data
  276. k, n = data[5], data[7]
  277. return data[9][:n-k-1]
  278. def get_residual(self):
  279. """Return weighted sum of squared residuals of the spline approximation.
  280. This is equivalent to::
  281. sum((w[i] * (y[i]-spl(x[i])))**2, axis=0)
  282. """
  283. return self._data[10]
  284. def integral(self, a, b):
  285. """ Return definite integral of the spline between two given points.
  286. Parameters
  287. ----------
  288. a : float
  289. Lower limit of integration.
  290. b : float
  291. Upper limit of integration.
  292. Returns
  293. -------
  294. integral : float
  295. The value of the definite integral of the spline between limits.
  296. Examples
  297. --------
  298. >>> from scipy.interpolate import UnivariateSpline
  299. >>> x = np.linspace(0, 3, 11)
  300. >>> y = x**2
  301. >>> spl = UnivariateSpline(x, y)
  302. >>> spl.integral(0, 3)
  303. 9.0
  304. which agrees with :math:`\\int x^2 dx = x^3 / 3` between the limits
  305. of 0 and 3.
  306. A caveat is that this routine assumes the spline to be zero outside of
  307. the data limits:
  308. >>> spl.integral(-1, 4)
  309. 9.0
  310. >>> spl.integral(-1, 0)
  311. 0.0
  312. """
  313. return dfitpack.splint(*(self._eval_args+(a, b)))
  314. def derivatives(self, x):
  315. """ Return all derivatives of the spline at the point x.
  316. Parameters
  317. ----------
  318. x : float
  319. The point to evaluate the derivatives at.
  320. Returns
  321. -------
  322. der : ndarray, shape(k+1,)
  323. Derivatives of the orders 0 to k.
  324. Examples
  325. --------
  326. >>> from scipy.interpolate import UnivariateSpline
  327. >>> x = np.linspace(0, 3, 11)
  328. >>> y = x**2
  329. >>> spl = UnivariateSpline(x, y)
  330. >>> spl.derivatives(1.5)
  331. array([2.25, 3.0, 2.0, 0])
  332. """
  333. d, ier = dfitpack.spalde(*(self._eval_args+(x,)))
  334. if not ier == 0:
  335. raise ValueError("Error code returned by spalde: %s" % ier)
  336. return d
  337. def roots(self):
  338. """ Return the zeros of the spline.
  339. Restriction: only cubic splines are supported by fitpack.
  340. """
  341. k = self._data[5]
  342. if k == 3:
  343. z, m, ier = dfitpack.sproot(*self._eval_args[:2])
  344. if not ier == 0:
  345. raise ValueError("Error code returned by spalde: %s" % ier)
  346. return z[:m]
  347. raise NotImplementedError('finding roots unsupported for '
  348. 'non-cubic splines')
  349. def derivative(self, n=1):
  350. """
  351. Construct a new spline representing the derivative of this spline.
  352. Parameters
  353. ----------
  354. n : int, optional
  355. Order of derivative to evaluate. Default: 1
  356. Returns
  357. -------
  358. spline : UnivariateSpline
  359. Spline of order k2=k-n representing the derivative of this
  360. spline.
  361. See Also
  362. --------
  363. splder, antiderivative
  364. Notes
  365. -----
  366. .. versionadded:: 0.13.0
  367. Examples
  368. --------
  369. This can be used for finding maxima of a curve:
  370. >>> from scipy.interpolate import UnivariateSpline
  371. >>> x = np.linspace(0, 10, 70)
  372. >>> y = np.sin(x)
  373. >>> spl = UnivariateSpline(x, y, k=4, s=0)
  374. Now, differentiate the spline and find the zeros of the
  375. derivative. (NB: `sproot` only works for order 3 splines, so we
  376. fit an order 4 spline):
  377. >>> spl.derivative().roots() / np.pi
  378. array([ 0.50000001, 1.5 , 2.49999998])
  379. This agrees well with roots :math:`\\pi/2 + n\\pi` of
  380. :math:`\\cos(x) = \\sin'(x)`.
  381. """
  382. tck = fitpack.splder(self._eval_args, n)
  383. return UnivariateSpline._from_tck(tck, self.ext)
  384. def antiderivative(self, n=1):
  385. """
  386. Construct a new spline representing the antiderivative of this spline.
  387. Parameters
  388. ----------
  389. n : int, optional
  390. Order of antiderivative to evaluate. Default: 1
  391. Returns
  392. -------
  393. spline : UnivariateSpline
  394. Spline of order k2=k+n representing the antiderivative of this
  395. spline.
  396. Notes
  397. -----
  398. .. versionadded:: 0.13.0
  399. See Also
  400. --------
  401. splantider, derivative
  402. Examples
  403. --------
  404. >>> from scipy.interpolate import UnivariateSpline
  405. >>> x = np.linspace(0, np.pi/2, 70)
  406. >>> y = 1 / np.sqrt(1 - 0.8*np.sin(x)**2)
  407. >>> spl = UnivariateSpline(x, y, s=0)
  408. The derivative is the inverse operation of the antiderivative,
  409. although some floating point error accumulates:
  410. >>> spl(1.7), spl.antiderivative().derivative()(1.7)
  411. (array(2.1565429877197317), array(2.1565429877201865))
  412. Antiderivative can be used to evaluate definite integrals:
  413. >>> ispl = spl.antiderivative()
  414. >>> ispl(np.pi/2) - ispl(0)
  415. 2.2572053588768486
  416. This is indeed an approximation to the complete elliptic integral
  417. :math:`K(m) = \\int_0^{\\pi/2} [1 - m\\sin^2 x]^{-1/2} dx`:
  418. >>> from scipy.special import ellipk
  419. >>> ellipk(0.8)
  420. 2.2572053268208538
  421. """
  422. tck = fitpack.splantider(self._eval_args, n)
  423. return UnivariateSpline._from_tck(tck, self.ext)
  424. class InterpolatedUnivariateSpline(UnivariateSpline):
  425. """
  426. One-dimensional interpolating spline for a given set of data points.
  427. Fits a spline y = spl(x) of degree `k` to the provided `x`, `y` data.
  428. Spline function passes through all provided points. Equivalent to
  429. `UnivariateSpline` with s=0.
  430. Parameters
  431. ----------
  432. x : (N,) array_like
  433. Input dimension of data points -- must be increasing
  434. y : (N,) array_like
  435. input dimension of data points
  436. w : (N,) array_like, optional
  437. Weights for spline fitting. Must be positive. If None (default),
  438. weights are all equal.
  439. bbox : (2,) array_like, optional
  440. 2-sequence specifying the boundary of the approximation interval. If
  441. None (default), ``bbox=[x[0], x[-1]]``.
  442. k : int, optional
  443. Degree of the smoothing spline. Must be 1 <= `k` <= 5.
  444. ext : int or str, optional
  445. Controls the extrapolation mode for elements
  446. not in the interval defined by the knot sequence.
  447. * if ext=0 or 'extrapolate', return the extrapolated value.
  448. * if ext=1 or 'zeros', return 0
  449. * if ext=2 or 'raise', raise a ValueError
  450. * if ext=3 of 'const', return the boundary value.
  451. The default value is 0.
  452. check_finite : bool, optional
  453. Whether to check that the input arrays contain only finite numbers.
  454. Disabling may give a performance gain, but may result in problems
  455. (crashes, non-termination or non-sensical results) if the inputs
  456. do contain infinities or NaNs.
  457. Default is False.
  458. See Also
  459. --------
  460. UnivariateSpline : Superclass -- allows knots to be selected by a
  461. smoothing condition
  462. LSQUnivariateSpline : spline for which knots are user-selected
  463. splrep : An older, non object-oriented wrapping of FITPACK
  464. splev, sproot, splint, spalde
  465. BivariateSpline : A similar class for two-dimensional spline interpolation
  466. Notes
  467. -----
  468. The number of data points must be larger than the spline degree `k`.
  469. Examples
  470. --------
  471. >>> import matplotlib.pyplot as plt
  472. >>> from scipy.interpolate import InterpolatedUnivariateSpline
  473. >>> x = np.linspace(-3, 3, 50)
  474. >>> y = np.exp(-x**2) + 0.1 * np.random.randn(50)
  475. >>> spl = InterpolatedUnivariateSpline(x, y)
  476. >>> plt.plot(x, y, 'ro', ms=5)
  477. >>> xs = np.linspace(-3, 3, 1000)
  478. >>> plt.plot(xs, spl(xs), 'g', lw=3, alpha=0.7)
  479. >>> plt.show()
  480. Notice that the ``spl(x)`` interpolates `y`:
  481. >>> spl.get_residual()
  482. 0.0
  483. """
  484. def __init__(self, x, y, w=None, bbox=[None]*2, k=3,
  485. ext=0, check_finite=False):
  486. if check_finite:
  487. w_finite = np.isfinite(w).all() if w is not None else True
  488. if (not np.isfinite(x).all() or not np.isfinite(y).all() or
  489. not w_finite):
  490. raise ValueError("Input must not contain NaNs or infs.")
  491. if not all(diff(x) > 0.0):
  492. raise ValueError('x must be strictly increasing')
  493. # _data == x,y,w,xb,xe,k,s,n,t,c,fp,fpint,nrdata,ier
  494. self._data = dfitpack.fpcurf0(x, y, k, w=w, xb=bbox[0],
  495. xe=bbox[1], s=0)
  496. self._reset_class()
  497. try:
  498. self.ext = _extrap_modes[ext]
  499. except KeyError:
  500. raise ValueError("Unknown extrapolation mode %s." % ext)
  501. _fpchec_error_string = """The input parameters have been rejected by fpchec. \
  502. This means that at least one of the following conditions is violated:
  503. 1) k+1 <= n-k-1 <= m
  504. 2) t(1) <= t(2) <= ... <= t(k+1)
  505. t(n-k) <= t(n-k+1) <= ... <= t(n)
  506. 3) t(k+1) < t(k+2) < ... < t(n-k)
  507. 4) t(k+1) <= x(i) <= t(n-k)
  508. 5) The conditions specified by Schoenberg and Whitney must hold
  509. for at least one subset of data points, i.e., there must be a
  510. subset of data points y(j) such that
  511. t(j) < y(j) < t(j+k+1), j=1,2,...,n-k-1
  512. """
  513. class LSQUnivariateSpline(UnivariateSpline):
  514. """
  515. One-dimensional spline with explicit internal knots.
  516. Fits a spline y = spl(x) of degree `k` to the provided `x`, `y` data. `t`
  517. specifies the internal knots of the spline
  518. Parameters
  519. ----------
  520. x : (N,) array_like
  521. Input dimension of data points -- must be increasing
  522. y : (N,) array_like
  523. Input dimension of data points
  524. t : (M,) array_like
  525. interior knots of the spline. Must be in ascending order and::
  526. bbox[0] < t[0] < ... < t[-1] < bbox[-1]
  527. w : (N,) array_like, optional
  528. weights for spline fitting. Must be positive. If None (default),
  529. weights are all equal.
  530. bbox : (2,) array_like, optional
  531. 2-sequence specifying the boundary of the approximation interval. If
  532. None (default), ``bbox = [x[0], x[-1]]``.
  533. k : int, optional
  534. Degree of the smoothing spline. Must be 1 <= `k` <= 5.
  535. Default is k=3, a cubic spline.
  536. ext : int or str, optional
  537. Controls the extrapolation mode for elements
  538. not in the interval defined by the knot sequence.
  539. * if ext=0 or 'extrapolate', return the extrapolated value.
  540. * if ext=1 or 'zeros', return 0
  541. * if ext=2 or 'raise', raise a ValueError
  542. * if ext=3 of 'const', return the boundary value.
  543. The default value is 0.
  544. check_finite : bool, optional
  545. Whether to check that the input arrays contain only finite numbers.
  546. Disabling may give a performance gain, but may result in problems
  547. (crashes, non-termination or non-sensical results) if the inputs
  548. do contain infinities or NaNs.
  549. Default is False.
  550. Raises
  551. ------
  552. ValueError
  553. If the interior knots do not satisfy the Schoenberg-Whitney conditions
  554. See Also
  555. --------
  556. UnivariateSpline : Superclass -- knots are specified by setting a
  557. smoothing condition
  558. InterpolatedUnivariateSpline : spline passing through all points
  559. splrep : An older, non object-oriented wrapping of FITPACK
  560. splev, sproot, splint, spalde
  561. BivariateSpline : A similar class for two-dimensional spline interpolation
  562. Notes
  563. -----
  564. The number of data points must be larger than the spline degree `k`.
  565. Knots `t` must satisfy the Schoenberg-Whitney conditions,
  566. i.e., there must be a subset of data points ``x[j]`` such that
  567. ``t[j] < x[j] < t[j+k+1]``, for ``j=0, 1,...,n-k-2``.
  568. Examples
  569. --------
  570. >>> from scipy.interpolate import LSQUnivariateSpline, UnivariateSpline
  571. >>> import matplotlib.pyplot as plt
  572. >>> x = np.linspace(-3, 3, 50)
  573. >>> y = np.exp(-x**2) + 0.1 * np.random.randn(50)
  574. Fit a smoothing spline with a pre-defined internal knots:
  575. >>> t = [-1, 0, 1]
  576. >>> spl = LSQUnivariateSpline(x, y, t)
  577. >>> xs = np.linspace(-3, 3, 1000)
  578. >>> plt.plot(x, y, 'ro', ms=5)
  579. >>> plt.plot(xs, spl(xs), 'g-', lw=3)
  580. >>> plt.show()
  581. Check the knot vector:
  582. >>> spl.get_knots()
  583. array([-3., -1., 0., 1., 3.])
  584. Constructing lsq spline using the knots from another spline:
  585. >>> x = np.arange(10)
  586. >>> s = UnivariateSpline(x, x, s=0)
  587. >>> s.get_knots()
  588. array([ 0., 2., 3., 4., 5., 6., 7., 9.])
  589. >>> knt = s.get_knots()
  590. >>> s1 = LSQUnivariateSpline(x, x, knt[1:-1]) # Chop 1st and last knot
  591. >>> s1.get_knots()
  592. array([ 0., 2., 3., 4., 5., 6., 7., 9.])
  593. """
  594. def __init__(self, x, y, t, w=None, bbox=[None]*2, k=3,
  595. ext=0, check_finite=False):
  596. if check_finite:
  597. w_finite = np.isfinite(w).all() if w is not None else True
  598. if (not np.isfinite(x).all() or not np.isfinite(y).all() or
  599. not w_finite or not np.isfinite(t).all()):
  600. raise ValueError("Input(s) must not contain NaNs or infs.")
  601. if not all(diff(x) > 0.0):
  602. raise ValueError('x must be strictly increasing')
  603. # _data == x,y,w,xb,xe,k,s,n,t,c,fp,fpint,nrdata,ier
  604. xb = bbox[0]
  605. xe = bbox[1]
  606. if xb is None:
  607. xb = x[0]
  608. if xe is None:
  609. xe = x[-1]
  610. t = concatenate(([xb]*(k+1), t, [xe]*(k+1)))
  611. n = len(t)
  612. if not alltrue(t[k+1:n-k]-t[k:n-k-1] > 0, axis=0):
  613. raise ValueError('Interior knots t must satisfy '
  614. 'Schoenberg-Whitney conditions')
  615. if not dfitpack.fpchec(x, t, k) == 0:
  616. raise ValueError(_fpchec_error_string)
  617. data = dfitpack.fpcurfm1(x, y, k, t, w=w, xb=xb, xe=xe)
  618. self._data = data[:-3] + (None, None, data[-1])
  619. self._reset_class()
  620. try:
  621. self.ext = _extrap_modes[ext]
  622. except KeyError:
  623. raise ValueError("Unknown extrapolation mode %s." % ext)
  624. # ############### Bivariate spline ####################
  625. class _BivariateSplineBase(object):
  626. """ Base class for Bivariate spline s(x,y) interpolation on the rectangle
  627. [xb,xe] x [yb, ye] calculated from a given set of data points
  628. (x,y,z).
  629. See Also
  630. --------
  631. bisplrep, bisplev : an older wrapping of FITPACK
  632. BivariateSpline :
  633. implementation of bivariate spline interpolation on a plane grid
  634. SphereBivariateSpline :
  635. implementation of bivariate spline interpolation on a spherical grid
  636. """
  637. def get_residual(self):
  638. """ Return weighted sum of squared residuals of the spline
  639. approximation: sum ((w[i]*(z[i]-s(x[i],y[i])))**2,axis=0)
  640. """
  641. return self.fp
  642. def get_knots(self):
  643. """ Return a tuple (tx,ty) where tx,ty contain knots positions
  644. of the spline with respect to x-, y-variable, respectively.
  645. The position of interior and additional knots are given as
  646. t[k+1:-k-1] and t[:k+1]=b, t[-k-1:]=e, respectively.
  647. """
  648. return self.tck[:2]
  649. def get_coeffs(self):
  650. """ Return spline coefficients."""
  651. return self.tck[2]
  652. def __call__(self, x, y, dx=0, dy=0, grid=True):
  653. """
  654. Evaluate the spline or its derivatives at given positions.
  655. Parameters
  656. ----------
  657. x, y : array_like
  658. Input coordinates.
  659. If `grid` is False, evaluate the spline at points ``(x[i],
  660. y[i]), i=0, ..., len(x)-1``. Standard Numpy broadcasting
  661. is obeyed.
  662. If `grid` is True: evaluate spline at the grid points
  663. defined by the coordinate arrays x, y. The arrays must be
  664. sorted to increasing order.
  665. Note that the axis ordering is inverted relative to
  666. the output of meshgrid.
  667. dx : int
  668. Order of x-derivative
  669. .. versionadded:: 0.14.0
  670. dy : int
  671. Order of y-derivative
  672. .. versionadded:: 0.14.0
  673. grid : bool
  674. Whether to evaluate the results on a grid spanned by the
  675. input arrays, or at points specified by the input arrays.
  676. .. versionadded:: 0.14.0
  677. """
  678. x = np.asarray(x)
  679. y = np.asarray(y)
  680. tx, ty, c = self.tck[:3]
  681. kx, ky = self.degrees
  682. if grid:
  683. if x.size == 0 or y.size == 0:
  684. return np.zeros((x.size, y.size), dtype=self.tck[2].dtype)
  685. if dx or dy:
  686. z, ier = dfitpack.parder(tx, ty, c, kx, ky, dx, dy, x, y)
  687. if not ier == 0:
  688. raise ValueError("Error code returned by parder: %s" % ier)
  689. else:
  690. z, ier = dfitpack.bispev(tx, ty, c, kx, ky, x, y)
  691. if not ier == 0:
  692. raise ValueError("Error code returned by bispev: %s" % ier)
  693. else:
  694. # standard Numpy broadcasting
  695. if x.shape != y.shape:
  696. x, y = np.broadcast_arrays(x, y)
  697. shape = x.shape
  698. x = x.ravel()
  699. y = y.ravel()
  700. if x.size == 0 or y.size == 0:
  701. return np.zeros(shape, dtype=self.tck[2].dtype)
  702. if dx or dy:
  703. z, ier = dfitpack.pardeu(tx, ty, c, kx, ky, dx, dy, x, y)
  704. if not ier == 0:
  705. raise ValueError("Error code returned by pardeu: %s" % ier)
  706. else:
  707. z, ier = dfitpack.bispeu(tx, ty, c, kx, ky, x, y)
  708. if not ier == 0:
  709. raise ValueError("Error code returned by bispeu: %s" % ier)
  710. z = z.reshape(shape)
  711. return z
  712. _surfit_messages = {1: """
  713. The required storage space exceeds the available storage space: nxest
  714. or nyest too small, or s too small.
  715. The weighted least-squares spline corresponds to the current set of
  716. knots.""",
  717. 2: """
  718. A theoretically impossible result was found during the iteration
  719. process for finding a smoothing spline with fp = s: s too small or
  720. badly chosen eps.
  721. Weighted sum of squared residuals does not satisfy abs(fp-s)/s < tol.""",
  722. 3: """
  723. the maximal number of iterations maxit (set to 20 by the program)
  724. allowed for finding a smoothing spline with fp=s has been reached:
  725. s too small.
  726. Weighted sum of squared residuals does not satisfy abs(fp-s)/s < tol.""",
  727. 4: """
  728. No more knots can be added because the number of b-spline coefficients
  729. (nx-kx-1)*(ny-ky-1) already exceeds the number of data points m:
  730. either s or m too small.
  731. The weighted least-squares spline corresponds to the current set of
  732. knots.""",
  733. 5: """
  734. No more knots can be added because the additional knot would (quasi)
  735. coincide with an old one: s too small or too large a weight to an
  736. inaccurate data point.
  737. The weighted least-squares spline corresponds to the current set of
  738. knots.""",
  739. 10: """
  740. Error on entry, no approximation returned. The following conditions
  741. must hold:
  742. xb<=x[i]<=xe, yb<=y[i]<=ye, w[i]>0, i=0..m-1
  743. If iopt==-1, then
  744. xb<tx[kx+1]<tx[kx+2]<...<tx[nx-kx-2]<xe
  745. yb<ty[ky+1]<ty[ky+2]<...<ty[ny-ky-2]<ye""",
  746. -3: """
  747. The coefficients of the spline returned have been computed as the
  748. minimal norm least-squares solution of a (numerically) rank deficient
  749. system (deficiency=%i). If deficiency is large, the results may be
  750. inaccurate. Deficiency may strongly depend on the value of eps."""
  751. }
  752. class BivariateSpline(_BivariateSplineBase):
  753. """
  754. Base class for bivariate splines.
  755. This describes a spline ``s(x, y)`` of degrees ``kx`` and ``ky`` on
  756. the rectangle ``[xb, xe] * [yb, ye]`` calculated from a given set
  757. of data points ``(x, y, z)``.
  758. This class is meant to be subclassed, not instantiated directly.
  759. To construct these splines, call either `SmoothBivariateSpline` or
  760. `LSQBivariateSpline`.
  761. See Also
  762. --------
  763. UnivariateSpline : a similar class for univariate spline interpolation
  764. SmoothBivariateSpline :
  765. to create a BivariateSpline through the given points
  766. LSQBivariateSpline :
  767. to create a BivariateSpline using weighted least-squares fitting
  768. SphereBivariateSpline :
  769. bivariate spline interpolation in spherical cooridinates
  770. bisplrep : older wrapping of FITPACK
  771. bisplev : older wrapping of FITPACK
  772. """
  773. @classmethod
  774. def _from_tck(cls, tck):
  775. """Construct a spline object from given tck and degree"""
  776. self = cls.__new__(cls)
  777. if len(tck) != 5:
  778. raise ValueError("tck should be a 5 element tuple of tx,"
  779. " ty, c, kx, ky")
  780. self.tck = tck[:3]
  781. self.degrees = tck[3:]
  782. return self
  783. def ev(self, xi, yi, dx=0, dy=0):
  784. """
  785. Evaluate the spline at points
  786. Returns the interpolated value at ``(xi[i], yi[i]),
  787. i=0,...,len(xi)-1``.
  788. Parameters
  789. ----------
  790. xi, yi : array_like
  791. Input coordinates. Standard Numpy broadcasting is obeyed.
  792. dx : int, optional
  793. Order of x-derivative
  794. .. versionadded:: 0.14.0
  795. dy : int, optional
  796. Order of y-derivative
  797. .. versionadded:: 0.14.0
  798. """
  799. return self.__call__(xi, yi, dx=dx, dy=dy, grid=False)
  800. def integral(self, xa, xb, ya, yb):
  801. """
  802. Evaluate the integral of the spline over area [xa,xb] x [ya,yb].
  803. Parameters
  804. ----------
  805. xa, xb : float
  806. The end-points of the x integration interval.
  807. ya, yb : float
  808. The end-points of the y integration interval.
  809. Returns
  810. -------
  811. integ : float
  812. The value of the resulting integral.
  813. """
  814. tx, ty, c = self.tck[:3]
  815. kx, ky = self.degrees
  816. return dfitpack.dblint(tx, ty, c, kx, ky, xa, xb, ya, yb)
  817. class SmoothBivariateSpline(BivariateSpline):
  818. """
  819. Smooth bivariate spline approximation.
  820. Parameters
  821. ----------
  822. x, y, z : array_like
  823. 1-D sequences of data points (order is not important).
  824. w : array_like, optional
  825. Positive 1-D sequence of weights, of same length as `x`, `y` and `z`.
  826. bbox : array_like, optional
  827. Sequence of length 4 specifying the boundary of the rectangular
  828. approximation domain. By default,
  829. ``bbox=[min(x,tx),max(x,tx), min(y,ty),max(y,ty)]``.
  830. kx, ky : ints, optional
  831. Degrees of the bivariate spline. Default is 3.
  832. s : float, optional
  833. Positive smoothing factor defined for estimation condition:
  834. ``sum((w[i]*(z[i]-s(x[i], y[i])))**2, axis=0) <= s``
  835. Default ``s=len(w)`` which should be a good value if ``1/w[i]`` is an
  836. estimate of the standard deviation of ``z[i]``.
  837. eps : float, optional
  838. A threshold for determining the effective rank of an over-determined
  839. linear system of equations. `eps` should have a value between 0 and 1,
  840. the default is 1e-16.
  841. See Also
  842. --------
  843. bisplrep : an older wrapping of FITPACK
  844. bisplev : an older wrapping of FITPACK
  845. UnivariateSpline : a similar class for univariate spline interpolation
  846. LSQUnivariateSpline : to create a BivariateSpline using weighted
  847. Notes
  848. -----
  849. The length of `x`, `y` and `z` should be at least ``(kx+1) * (ky+1)``.
  850. """
  851. def __init__(self, x, y, z, w=None, bbox=[None] * 4, kx=3, ky=3, s=None,
  852. eps=None):
  853. xb, xe, yb, ye = bbox
  854. nx, tx, ny, ty, c, fp, wrk1, ier = dfitpack.surfit_smth(x, y, z, w,
  855. xb, xe, yb,
  856. ye, kx, ky,
  857. s=s, eps=eps,
  858. lwrk2=1)
  859. if ier > 10: # lwrk2 was to small, re-run
  860. nx, tx, ny, ty, c, fp, wrk1, ier = dfitpack.surfit_smth(x, y, z, w,
  861. xb, xe, yb,
  862. ye, kx, ky,
  863. s=s,
  864. eps=eps,
  865. lwrk2=ier)
  866. if ier in [0, -1, -2]: # normal return
  867. pass
  868. else:
  869. message = _surfit_messages.get(ier, 'ier=%s' % (ier))
  870. warnings.warn(message)
  871. self.fp = fp
  872. self.tck = tx[:nx], ty[:ny], c[:(nx-kx-1)*(ny-ky-1)]
  873. self.degrees = kx, ky
  874. class LSQBivariateSpline(BivariateSpline):
  875. """
  876. Weighted least-squares bivariate spline approximation.
  877. Parameters
  878. ----------
  879. x, y, z : array_like
  880. 1-D sequences of data points (order is not important).
  881. tx, ty : array_like
  882. Strictly ordered 1-D sequences of knots coordinates.
  883. w : array_like, optional
  884. Positive 1-D array of weights, of the same length as `x`, `y` and `z`.
  885. bbox : (4,) array_like, optional
  886. Sequence of length 4 specifying the boundary of the rectangular
  887. approximation domain. By default,
  888. ``bbox=[min(x,tx),max(x,tx), min(y,ty),max(y,ty)]``.
  889. kx, ky : ints, optional
  890. Degrees of the bivariate spline. Default is 3.
  891. eps : float, optional
  892. A threshold for determining the effective rank of an over-determined
  893. linear system of equations. `eps` should have a value between 0 and 1,
  894. the default is 1e-16.
  895. See Also
  896. --------
  897. bisplrep : an older wrapping of FITPACK
  898. bisplev : an older wrapping of FITPACK
  899. UnivariateSpline : a similar class for univariate spline interpolation
  900. SmoothBivariateSpline : create a smoothing BivariateSpline
  901. Notes
  902. -----
  903. The length of `x`, `y` and `z` should be at least ``(kx+1) * (ky+1)``.
  904. """
  905. def __init__(self, x, y, z, tx, ty, w=None, bbox=[None]*4, kx=3, ky=3,
  906. eps=None):
  907. nx = 2*kx+2+len(tx)
  908. ny = 2*ky+2+len(ty)
  909. tx1 = zeros((nx,), float)
  910. ty1 = zeros((ny,), float)
  911. tx1[kx+1:nx-kx-1] = tx
  912. ty1[ky+1:ny-ky-1] = ty
  913. xb, xe, yb, ye = bbox
  914. tx1, ty1, c, fp, ier = dfitpack.surfit_lsq(x, y, z, tx1, ty1, w,
  915. xb, xe, yb, ye,
  916. kx, ky, eps, lwrk2=1)
  917. if ier > 10:
  918. tx1, ty1, c, fp, ier = dfitpack.surfit_lsq(x, y, z, tx1, ty1, w,
  919. xb, xe, yb, ye,
  920. kx, ky, eps, lwrk2=ier)
  921. if ier in [0, -1, -2]: # normal return
  922. pass
  923. else:
  924. if ier < -2:
  925. deficiency = (nx-kx-1)*(ny-ky-1)+ier
  926. message = _surfit_messages.get(-3) % (deficiency)
  927. else:
  928. message = _surfit_messages.get(ier, 'ier=%s' % (ier))
  929. warnings.warn(message)
  930. self.fp = fp
  931. self.tck = tx1, ty1, c
  932. self.degrees = kx, ky
  933. class RectBivariateSpline(BivariateSpline):
  934. """
  935. Bivariate spline approximation over a rectangular mesh.
  936. Can be used for both smoothing and interpolating data.
  937. Parameters
  938. ----------
  939. x,y : array_like
  940. 1-D arrays of coordinates in strictly ascending order.
  941. z : array_like
  942. 2-D array of data with shape (x.size,y.size).
  943. bbox : array_like, optional
  944. Sequence of length 4 specifying the boundary of the rectangular
  945. approximation domain. By default,
  946. ``bbox=[min(x,tx),max(x,tx), min(y,ty),max(y,ty)]``.
  947. kx, ky : ints, optional
  948. Degrees of the bivariate spline. Default is 3.
  949. s : float, optional
  950. Positive smoothing factor defined for estimation condition:
  951. ``sum((w[i]*(z[i]-s(x[i], y[i])))**2, axis=0) <= s``
  952. Default is ``s=0``, which is for interpolation.
  953. See Also
  954. --------
  955. SmoothBivariateSpline : a smoothing bivariate spline for scattered data
  956. bisplrep : an older wrapping of FITPACK
  957. bisplev : an older wrapping of FITPACK
  958. UnivariateSpline : a similar class for univariate spline interpolation
  959. """
  960. def __init__(self, x, y, z, bbox=[None] * 4, kx=3, ky=3, s=0):
  961. x, y = ravel(x), ravel(y)
  962. if not all(diff(x) > 0.0):
  963. raise ValueError('x must be strictly increasing')
  964. if not all(diff(y) > 0.0):
  965. raise ValueError('y must be strictly increasing')
  966. if not ((x.min() == x[0]) and (x.max() == x[-1])):
  967. raise ValueError('x must be strictly ascending')
  968. if not ((y.min() == y[0]) and (y.max() == y[-1])):
  969. raise ValueError('y must be strictly ascending')
  970. if not x.size == z.shape[0]:
  971. raise ValueError('x dimension of z must have same number of '
  972. 'elements as x')
  973. if not y.size == z.shape[1]:
  974. raise ValueError('y dimension of z must have same number of '
  975. 'elements as y')
  976. z = ravel(z)
  977. xb, xe, yb, ye = bbox
  978. nx, tx, ny, ty, c, fp, ier = dfitpack.regrid_smth(x, y, z, xb, xe, yb,
  979. ye, kx, ky, s)
  980. if ier not in [0, -1, -2]:
  981. msg = _surfit_messages.get(ier, 'ier=%s' % (ier))
  982. raise ValueError(msg)
  983. self.fp = fp
  984. self.tck = tx[:nx], ty[:ny], c[:(nx - kx - 1) * (ny - ky - 1)]
  985. self.degrees = kx, ky
  986. _spherefit_messages = _surfit_messages.copy()
  987. _spherefit_messages[10] = """
  988. ERROR. On entry, the input data are controlled on validity. The following
  989. restrictions must be satisfied:
  990. -1<=iopt<=1, m>=2, ntest>=8 ,npest >=8, 0<eps<1,
  991. 0<=teta(i)<=pi, 0<=phi(i)<=2*pi, w(i)>0, i=1,...,m
  992. lwrk1 >= 185+52*v+10*u+14*u*v+8*(u-1)*v**2+8*m
  993. kwrk >= m+(ntest-7)*(npest-7)
  994. if iopt=-1: 8<=nt<=ntest , 9<=np<=npest
  995. 0<tt(5)<tt(6)<...<tt(nt-4)<pi
  996. 0<tp(5)<tp(6)<...<tp(np-4)<2*pi
  997. if iopt>=0: s>=0
  998. if one of these conditions is found to be violated,control
  999. is immediately repassed to the calling program. in that
  1000. case there is no approximation returned."""
  1001. _spherefit_messages[-3] = """
  1002. WARNING. The coefficients of the spline returned have been computed as the
  1003. minimal norm least-squares solution of a (numerically) rank
  1004. deficient system (deficiency=%i, rank=%i). Especially if the rank
  1005. deficiency, which is computed by 6+(nt-8)*(np-7)+ier, is large,
  1006. the results may be inaccurate. They could also seriously depend on
  1007. the value of eps."""
  1008. class SphereBivariateSpline(_BivariateSplineBase):
  1009. """
  1010. Bivariate spline s(x,y) of degrees 3 on a sphere, calculated from a
  1011. given set of data points (theta,phi,r).
  1012. .. versionadded:: 0.11.0
  1013. See Also
  1014. --------
  1015. bisplrep, bisplev : an older wrapping of FITPACK
  1016. UnivariateSpline : a similar class for univariate spline interpolation
  1017. SmoothUnivariateSpline :
  1018. to create a BivariateSpline through the given points
  1019. LSQUnivariateSpline :
  1020. to create a BivariateSpline using weighted least-squares fitting
  1021. """
  1022. def __call__(self, theta, phi, dtheta=0, dphi=0, grid=True):
  1023. """
  1024. Evaluate the spline or its derivatives at given positions.
  1025. Parameters
  1026. ----------
  1027. theta, phi : array_like
  1028. Input coordinates.
  1029. If `grid` is False, evaluate the spline at points
  1030. ``(theta[i], phi[i]), i=0, ..., len(x)-1``. Standard
  1031. Numpy broadcasting is obeyed.
  1032. If `grid` is True: evaluate spline at the grid points
  1033. defined by the coordinate arrays theta, phi. The arrays
  1034. must be sorted to increasing order.
  1035. dtheta : int, optional
  1036. Order of theta-derivative
  1037. .. versionadded:: 0.14.0
  1038. dphi : int
  1039. Order of phi-derivative
  1040. .. versionadded:: 0.14.0
  1041. grid : bool
  1042. Whether to evaluate the results on a grid spanned by the
  1043. input arrays, or at points specified by the input arrays.
  1044. .. versionadded:: 0.14.0
  1045. """
  1046. theta = np.asarray(theta)
  1047. phi = np.asarray(phi)
  1048. if theta.size > 0 and (theta.min() < 0. or theta.max() > np.pi):
  1049. raise ValueError("requested theta out of bounds.")
  1050. if phi.size > 0 and (phi.min() < 0. or phi.max() > 2. * np.pi):
  1051. raise ValueError("requested phi out of bounds.")
  1052. return _BivariateSplineBase.__call__(self, theta, phi,
  1053. dx=dtheta, dy=dphi, grid=grid)
  1054. def ev(self, theta, phi, dtheta=0, dphi=0):
  1055. """
  1056. Evaluate the spline at points
  1057. Returns the interpolated value at ``(theta[i], phi[i]),
  1058. i=0,...,len(theta)-1``.
  1059. Parameters
  1060. ----------
  1061. theta, phi : array_like
  1062. Input coordinates. Standard Numpy broadcasting is obeyed.
  1063. dtheta : int, optional
  1064. Order of theta-derivative
  1065. .. versionadded:: 0.14.0
  1066. dphi : int, optional
  1067. Order of phi-derivative
  1068. .. versionadded:: 0.14.0
  1069. """
  1070. return self.__call__(theta, phi, dtheta=dtheta, dphi=dphi, grid=False)
  1071. class SmoothSphereBivariateSpline(SphereBivariateSpline):
  1072. """
  1073. Smooth bivariate spline approximation in spherical coordinates.
  1074. .. versionadded:: 0.11.0
  1075. Parameters
  1076. ----------
  1077. theta, phi, r : array_like
  1078. 1-D sequences of data points (order is not important). Coordinates
  1079. must be given in radians. Theta must lie within the interval (0, pi),
  1080. and phi must lie within the interval (0, 2pi).
  1081. w : array_like, optional
  1082. Positive 1-D sequence of weights.
  1083. s : float, optional
  1084. Positive smoothing factor defined for estimation condition:
  1085. ``sum((w(i)*(r(i) - s(theta(i), phi(i))))**2, axis=0) <= s``
  1086. Default ``s=len(w)`` which should be a good value if 1/w[i] is an
  1087. estimate of the standard deviation of r[i].
  1088. eps : float, optional
  1089. A threshold for determining the effective rank of an over-determined
  1090. linear system of equations. `eps` should have a value between 0 and 1,
  1091. the default is 1e-16.
  1092. Notes
  1093. -----
  1094. For more information, see the FITPACK_ site about this function.
  1095. .. _FITPACK: http://www.netlib.org/dierckx/sphere.f
  1096. Examples
  1097. --------
  1098. Suppose we have global data on a coarse grid (the input data does not
  1099. have to be on a grid):
  1100. >>> theta = np.linspace(0., np.pi, 7)
  1101. >>> phi = np.linspace(0., 2*np.pi, 9)
  1102. >>> data = np.empty((theta.shape[0], phi.shape[0]))
  1103. >>> data[:,0], data[0,:], data[-1,:] = 0., 0., 0.
  1104. >>> data[1:-1,1], data[1:-1,-1] = 1., 1.
  1105. >>> data[1,1:-1], data[-2,1:-1] = 1., 1.
  1106. >>> data[2:-2,2], data[2:-2,-2] = 2., 2.
  1107. >>> data[2,2:-2], data[-3,2:-2] = 2., 2.
  1108. >>> data[3,3:-2] = 3.
  1109. >>> data = np.roll(data, 4, 1)
  1110. We need to set up the interpolator object
  1111. >>> lats, lons = np.meshgrid(theta, phi)
  1112. >>> from scipy.interpolate import SmoothSphereBivariateSpline
  1113. >>> lut = SmoothSphereBivariateSpline(lats.ravel(), lons.ravel(),
  1114. ... data.T.ravel(), s=3.5)
  1115. As a first test, we'll see what the algorithm returns when run on the
  1116. input coordinates
  1117. >>> data_orig = lut(theta, phi)
  1118. Finally we interpolate the data to a finer grid
  1119. >>> fine_lats = np.linspace(0., np.pi, 70)
  1120. >>> fine_lons = np.linspace(0., 2 * np.pi, 90)
  1121. >>> data_smth = lut(fine_lats, fine_lons)
  1122. >>> import matplotlib.pyplot as plt
  1123. >>> fig = plt.figure()
  1124. >>> ax1 = fig.add_subplot(131)
  1125. >>> ax1.imshow(data, interpolation='nearest')
  1126. >>> ax2 = fig.add_subplot(132)
  1127. >>> ax2.imshow(data_orig, interpolation='nearest')
  1128. >>> ax3 = fig.add_subplot(133)
  1129. >>> ax3.imshow(data_smth, interpolation='nearest')
  1130. >>> plt.show()
  1131. """
  1132. def __init__(self, theta, phi, r, w=None, s=0., eps=1E-16):
  1133. if np.issubclass_(w, float):
  1134. w = ones(len(theta)) * w
  1135. nt_, tt_, np_, tp_, c, fp, ier = dfitpack.spherfit_smth(theta, phi,
  1136. r, w=w, s=s,
  1137. eps=eps)
  1138. if ier not in [0, -1, -2]:
  1139. message = _spherefit_messages.get(ier, 'ier=%s' % (ier))
  1140. raise ValueError(message)
  1141. self.fp = fp
  1142. self.tck = tt_[:nt_], tp_[:np_], c[:(nt_ - 4) * (np_ - 4)]
  1143. self.degrees = (3, 3)
  1144. class LSQSphereBivariateSpline(SphereBivariateSpline):
  1145. """
  1146. Weighted least-squares bivariate spline approximation in spherical
  1147. coordinates.
  1148. .. versionadded:: 0.11.0
  1149. Parameters
  1150. ----------
  1151. theta, phi, r : array_like
  1152. 1-D sequences of data points (order is not important). Coordinates
  1153. must be given in radians. Theta must lie within the interval (0, pi),
  1154. and phi must lie within the interval (0, 2pi).
  1155. tt, tp : array_like
  1156. Strictly ordered 1-D sequences of knots coordinates.
  1157. Coordinates must satisfy ``0 < tt[i] < pi``, ``0 < tp[i] < 2*pi``.
  1158. w : array_like, optional
  1159. Positive 1-D sequence of weights, of the same length as `theta`, `phi`
  1160. and `r`.
  1161. eps : float, optional
  1162. A threshold for determining the effective rank of an over-determined
  1163. linear system of equations. `eps` should have a value between 0 and 1,
  1164. the default is 1e-16.
  1165. Notes
  1166. -----
  1167. For more information, see the FITPACK_ site about this function.
  1168. .. _FITPACK: http://www.netlib.org/dierckx/sphere.f
  1169. Examples
  1170. --------
  1171. Suppose we have global data on a coarse grid (the input data does not
  1172. have to be on a grid):
  1173. >>> theta = np.linspace(0., np.pi, 7)
  1174. >>> phi = np.linspace(0., 2*np.pi, 9)
  1175. >>> data = np.empty((theta.shape[0], phi.shape[0]))
  1176. >>> data[:,0], data[0,:], data[-1,:] = 0., 0., 0.
  1177. >>> data[1:-1,1], data[1:-1,-1] = 1., 1.
  1178. >>> data[1,1:-1], data[-2,1:-1] = 1., 1.
  1179. >>> data[2:-2,2], data[2:-2,-2] = 2., 2.
  1180. >>> data[2,2:-2], data[-3,2:-2] = 2., 2.
  1181. >>> data[3,3:-2] = 3.
  1182. >>> data = np.roll(data, 4, 1)
  1183. We need to set up the interpolator object. Here, we must also specify the
  1184. coordinates of the knots to use.
  1185. >>> lats, lons = np.meshgrid(theta, phi)
  1186. >>> knotst, knotsp = theta.copy(), phi.copy()
  1187. >>> knotst[0] += .0001
  1188. >>> knotst[-1] -= .0001
  1189. >>> knotsp[0] += .0001
  1190. >>> knotsp[-1] -= .0001
  1191. >>> from scipy.interpolate import LSQSphereBivariateSpline
  1192. >>> lut = LSQSphereBivariateSpline(lats.ravel(), lons.ravel(),
  1193. ... data.T.ravel(), knotst, knotsp)
  1194. As a first test, we'll see what the algorithm returns when run on the
  1195. input coordinates
  1196. >>> data_orig = lut(theta, phi)
  1197. Finally we interpolate the data to a finer grid
  1198. >>> fine_lats = np.linspace(0., np.pi, 70)
  1199. >>> fine_lons = np.linspace(0., 2*np.pi, 90)
  1200. >>> data_lsq = lut(fine_lats, fine_lons)
  1201. >>> import matplotlib.pyplot as plt
  1202. >>> fig = plt.figure()
  1203. >>> ax1 = fig.add_subplot(131)
  1204. >>> ax1.imshow(data, interpolation='nearest')
  1205. >>> ax2 = fig.add_subplot(132)
  1206. >>> ax2.imshow(data_orig, interpolation='nearest')
  1207. >>> ax3 = fig.add_subplot(133)
  1208. >>> ax3.imshow(data_lsq, interpolation='nearest')
  1209. >>> plt.show()
  1210. """
  1211. def __init__(self, theta, phi, r, tt, tp, w=None, eps=1E-16):
  1212. if np.issubclass_(w, float):
  1213. w = ones(len(theta)) * w
  1214. nt_, np_ = 8 + len(tt), 8 + len(tp)
  1215. tt_, tp_ = zeros((nt_,), float), zeros((np_,), float)
  1216. tt_[4:-4], tp_[4:-4] = tt, tp
  1217. tt_[-4:], tp_[-4:] = np.pi, 2. * np.pi
  1218. tt_, tp_, c, fp, ier = dfitpack.spherfit_lsq(theta, phi, r, tt_, tp_,
  1219. w=w, eps=eps)
  1220. if ier < -2:
  1221. deficiency = 6 + (nt_ - 8) * (np_ - 7) + ier
  1222. message = _spherefit_messages.get(-3) % (deficiency, -ier)
  1223. warnings.warn(message, stacklevel=3)
  1224. elif ier not in [0, -1, -2]:
  1225. message = _spherefit_messages.get(ier, 'ier=%s' % (ier))
  1226. raise ValueError(message)
  1227. self.fp = fp
  1228. self.tck = tt_, tp_, c
  1229. self.degrees = (3, 3)
  1230. _spfit_messages = _surfit_messages.copy()
  1231. _spfit_messages[10] = """
  1232. ERROR: on entry, the input data are controlled on validity
  1233. the following restrictions must be satisfied.
  1234. -1<=iopt(1)<=1, 0<=iopt(2)<=1, 0<=iopt(3)<=1,
  1235. -1<=ider(1)<=1, 0<=ider(2)<=1, ider(2)=0 if iopt(2)=0.
  1236. -1<=ider(3)<=1, 0<=ider(4)<=1, ider(4)=0 if iopt(3)=0.
  1237. mu >= mumin (see above), mv >= 4, nuest >=8, nvest >= 8,
  1238. kwrk>=5+mu+mv+nuest+nvest,
  1239. lwrk >= 12+nuest*(mv+nvest+3)+nvest*24+4*mu+8*mv+max(nuest,mv+nvest)
  1240. 0< u(i-1)<u(i)< pi,i=2,..,mu,
  1241. -pi<=v(1)< pi, v(1)<v(i-1)<v(i)<v(1)+2*pi, i=3,...,mv
  1242. if iopt(1)=-1: 8<=nu<=min(nuest,mu+6+iopt(2)+iopt(3))
  1243. 0<tu(5)<tu(6)<...<tu(nu-4)< pi
  1244. 8<=nv<=min(nvest,mv+7)
  1245. v(1)<tv(5)<tv(6)<...<tv(nv-4)<v(1)+2*pi
  1246. the schoenberg-whitney conditions, i.e. there must be
  1247. subset of grid co-ordinates uu(p) and vv(q) such that
  1248. tu(p) < uu(p) < tu(p+4) ,p=1,...,nu-4
  1249. (iopt(2)=1 and iopt(3)=1 also count for a uu-value
  1250. tv(q) < vv(q) < tv(q+4) ,q=1,...,nv-4
  1251. (vv(q) is either a value v(j) or v(j)+2*pi)
  1252. if iopt(1)>=0: s>=0
  1253. if s=0: nuest>=mu+6+iopt(2)+iopt(3), nvest>=mv+7
  1254. if one of these conditions is found to be violated,control is
  1255. immediately repassed to the calling program. in that case there is no
  1256. approximation returned."""
  1257. class RectSphereBivariateSpline(SphereBivariateSpline):
  1258. """
  1259. Bivariate spline approximation over a rectangular mesh on a sphere.
  1260. Can be used for smoothing data.
  1261. .. versionadded:: 0.11.0
  1262. Parameters
  1263. ----------
  1264. u : array_like
  1265. 1-D array of latitude coordinates in strictly ascending order.
  1266. Coordinates must be given in radians and lie within the interval
  1267. (0, pi).
  1268. v : array_like
  1269. 1-D array of longitude coordinates in strictly ascending order.
  1270. Coordinates must be given in radians. First element (v[0]) must lie
  1271. within the interval [-pi, pi). Last element (v[-1]) must satisfy
  1272. v[-1] <= v[0] + 2*pi.
  1273. r : array_like
  1274. 2-D array of data with shape ``(u.size, v.size)``.
  1275. s : float, optional
  1276. Positive smoothing factor defined for estimation condition
  1277. (``s=0`` is for interpolation).
  1278. pole_continuity : bool or (bool, bool), optional
  1279. Order of continuity at the poles ``u=0`` (``pole_continuity[0]``) and
  1280. ``u=pi`` (``pole_continuity[1]``). The order of continuity at the pole
  1281. will be 1 or 0 when this is True or False, respectively.
  1282. Defaults to False.
  1283. pole_values : float or (float, float), optional
  1284. Data values at the poles ``u=0`` and ``u=pi``. Either the whole
  1285. parameter or each individual element can be None. Defaults to None.
  1286. pole_exact : bool or (bool, bool), optional
  1287. Data value exactness at the poles ``u=0`` and ``u=pi``. If True, the
  1288. value is considered to be the right function value, and it will be
  1289. fitted exactly. If False, the value will be considered to be a data
  1290. value just like the other data values. Defaults to False.
  1291. pole_flat : bool or (bool, bool), optional
  1292. For the poles at ``u=0`` and ``u=pi``, specify whether or not the
  1293. approximation has vanishing derivatives. Defaults to False.
  1294. See Also
  1295. --------
  1296. RectBivariateSpline : bivariate spline approximation over a rectangular
  1297. mesh
  1298. Notes
  1299. -----
  1300. Currently, only the smoothing spline approximation (``iopt[0] = 0`` and
  1301. ``iopt[0] = 1`` in the FITPACK routine) is supported. The exact
  1302. least-squares spline approximation is not implemented yet.
  1303. When actually performing the interpolation, the requested `v` values must
  1304. lie within the same length 2pi interval that the original `v` values were
  1305. chosen from.
  1306. For more information, see the FITPACK_ site about this function.
  1307. .. _FITPACK: http://www.netlib.org/dierckx/spgrid.f
  1308. Examples
  1309. --------
  1310. Suppose we have global data on a coarse grid
  1311. >>> lats = np.linspace(10, 170, 9) * np.pi / 180.
  1312. >>> lons = np.linspace(0, 350, 18) * np.pi / 180.
  1313. >>> data = np.dot(np.atleast_2d(90. - np.linspace(-80., 80., 18)).T,
  1314. ... np.atleast_2d(180. - np.abs(np.linspace(0., 350., 9)))).T
  1315. We want to interpolate it to a global one-degree grid
  1316. >>> new_lats = np.linspace(1, 180, 180) * np.pi / 180
  1317. >>> new_lons = np.linspace(1, 360, 360) * np.pi / 180
  1318. >>> new_lats, new_lons = np.meshgrid(new_lats, new_lons)
  1319. We need to set up the interpolator object
  1320. >>> from scipy.interpolate import RectSphereBivariateSpline
  1321. >>> lut = RectSphereBivariateSpline(lats, lons, data)
  1322. Finally we interpolate the data. The `RectSphereBivariateSpline` object
  1323. only takes 1-D arrays as input, therefore we need to do some reshaping.
  1324. >>> data_interp = lut.ev(new_lats.ravel(),
  1325. ... new_lons.ravel()).reshape((360, 180)).T
  1326. Looking at the original and the interpolated data, one can see that the
  1327. interpolant reproduces the original data very well:
  1328. >>> import matplotlib.pyplot as plt
  1329. >>> fig = plt.figure()
  1330. >>> ax1 = fig.add_subplot(211)
  1331. >>> ax1.imshow(data, interpolation='nearest')
  1332. >>> ax2 = fig.add_subplot(212)
  1333. >>> ax2.imshow(data_interp, interpolation='nearest')
  1334. >>> plt.show()
  1335. Choosing the optimal value of ``s`` can be a delicate task. Recommended
  1336. values for ``s`` depend on the accuracy of the data values. If the user
  1337. has an idea of the statistical errors on the data, she can also find a
  1338. proper estimate for ``s``. By assuming that, if she specifies the
  1339. right ``s``, the interpolator will use a spline ``f(u,v)`` which exactly
  1340. reproduces the function underlying the data, she can evaluate
  1341. ``sum((r(i,j)-s(u(i),v(j)))**2)`` to find a good estimate for this ``s``.
  1342. For example, if she knows that the statistical errors on her
  1343. ``r(i,j)``-values are not greater than 0.1, she may expect that a good
  1344. ``s`` should have a value not larger than ``u.size * v.size * (0.1)**2``.
  1345. If nothing is known about the statistical error in ``r(i,j)``, ``s`` must
  1346. be determined by trial and error. The best is then to start with a very
  1347. large value of ``s`` (to determine the least-squares polynomial and the
  1348. corresponding upper bound ``fp0`` for ``s``) and then to progressively
  1349. decrease the value of ``s`` (say by a factor 10 in the beginning, i.e.
  1350. ``s = fp0 / 10, fp0 / 100, ...`` and more carefully as the approximation
  1351. shows more detail) to obtain closer fits.
  1352. The interpolation results for different values of ``s`` give some insight
  1353. into this process:
  1354. >>> fig2 = plt.figure()
  1355. >>> s = [3e9, 2e9, 1e9, 1e8]
  1356. >>> for ii in range(len(s)):
  1357. ... lut = RectSphereBivariateSpline(lats, lons, data, s=s[ii])
  1358. ... data_interp = lut.ev(new_lats.ravel(),
  1359. ... new_lons.ravel()).reshape((360, 180)).T
  1360. ... ax = fig2.add_subplot(2, 2, ii+1)
  1361. ... ax.imshow(data_interp, interpolation='nearest')
  1362. ... ax.set_title("s = %g" % s[ii])
  1363. >>> plt.show()
  1364. """
  1365. def __init__(self, u, v, r, s=0., pole_continuity=False, pole_values=None,
  1366. pole_exact=False, pole_flat=False):
  1367. iopt = np.array([0, 0, 0], dtype=int)
  1368. ider = np.array([-1, 0, -1, 0], dtype=int)
  1369. if pole_values is None:
  1370. pole_values = (None, None)
  1371. elif isinstance(pole_values, (float, np.float32, np.float64)):
  1372. pole_values = (pole_values, pole_values)
  1373. if isinstance(pole_continuity, bool):
  1374. pole_continuity = (pole_continuity, pole_continuity)
  1375. if isinstance(pole_exact, bool):
  1376. pole_exact = (pole_exact, pole_exact)
  1377. if isinstance(pole_flat, bool):
  1378. pole_flat = (pole_flat, pole_flat)
  1379. r0, r1 = pole_values
  1380. iopt[1:] = pole_continuity
  1381. if r0 is None:
  1382. ider[0] = -1
  1383. else:
  1384. ider[0] = pole_exact[0]
  1385. if r1 is None:
  1386. ider[2] = -1
  1387. else:
  1388. ider[2] = pole_exact[1]
  1389. ider[1], ider[3] = pole_flat
  1390. u, v = np.ravel(u), np.ravel(v)
  1391. if not np.all(np.diff(u) > 0.0):
  1392. raise ValueError('u must be strictly increasing')
  1393. if not np.all(np.diff(v) > 0.0):
  1394. raise ValueError('v must be strictly increasing')
  1395. if not u.size == r.shape[0]:
  1396. raise ValueError('u dimension of r must have same number of '
  1397. 'elements as u')
  1398. if not v.size == r.shape[1]:
  1399. raise ValueError('v dimension of r must have same number of '
  1400. 'elements as v')
  1401. if pole_continuity[1] is False and pole_flat[1] is True:
  1402. raise ValueError('if pole_continuity is False, so must be '
  1403. 'pole_flat')
  1404. if pole_continuity[0] is False and pole_flat[0] is True:
  1405. raise ValueError('if pole_continuity is False, so must be '
  1406. 'pole_flat')
  1407. r = np.ravel(r)
  1408. nu, tu, nv, tv, c, fp, ier = dfitpack.regrid_smth_spher(iopt, ider,
  1409. u.copy(), v.copy(), r.copy(), r0, r1, s)
  1410. if ier not in [0, -1, -2]:
  1411. msg = _spfit_messages.get(ier, 'ier=%s' % (ier))
  1412. raise ValueError(msg)
  1413. self.fp = fp
  1414. self.tck = tu[:nu], tv[:nv], c[:(nu - 4) * (nv-4)]
  1415. self.degrees = (3, 3)