1"""
2fitpack (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
6See
7 https://web.archive.org/web/20010524124604/http://www.cs.kuleuven.ac.be:80/cwis/research/nalag/research/topics/fitpack.html
8or
9 http://www.netlib.org/dierckx/
10
11Copyright 2002 Pearu Peterson all rights reserved,
12Pearu Peterson <pearu@cens.ioc.ee>
13Permission to use, modify, and distribute this software is given under the
14terms of the SciPy (BSD style) license. See LICENSE.txt that came with
15this distribution for specifics.
16
17NO WARRANTY IS EXPRESSED OR IMPLIED.  USE AT YOUR OWN RISK.
18
19TODO: Make interfaces to the following fitpack functions:
20    For univariate splines: cocosp, concon, fourco, insert
21    For bivariate splines: profil, regrid, parsur, surev
22"""
23
24__all__ = ['splrep', 'splprep', 'splev', 'splint', 'sproot', 'spalde',
25           'bisplrep', 'bisplev', 'insert', 'splder', 'splantider']
26
27import warnings
28import numpy as np
29from . import _fitpack
30from numpy import (atleast_1d, array, ones, zeros, sqrt, ravel, transpose,
31                   empty, iinfo, asarray)
32
33# Try to replace _fitpack interface with
34#  f2py-generated version
35from . import dfitpack
36
37
38dfitpack_int = dfitpack.types.intvar.dtype
39
40
41def _int_overflow(x, msg=None):
42    """Cast the value to an dfitpack_int and raise an OverflowError if the value
43    cannot fit.
44    """
45    if x > iinfo(dfitpack_int).max:
46        if msg is None:
47            msg = '%r cannot fit into an %r' % (x, dfitpack_int)
48        raise OverflowError(msg)
49    return dfitpack_int.type(x)
50
51
52_iermess = {
53    0: ["The spline has a residual sum of squares fp such that "
54        "abs(fp-s)/s<=0.001", None],
55    -1: ["The spline is an interpolating spline (fp=0)", None],
56    -2: ["The spline is weighted least-squares polynomial of degree k.\n"
57         "fp gives the upper bound fp0 for the smoothing factor s", None],
58    1: ["The required storage space exceeds the available storage space.\n"
59        "Probable causes: data (x,y) size is too small or smoothing parameter"
60        "\ns is too small (fp>s).", ValueError],
61    2: ["A theoretically impossible result when finding a smoothing spline\n"
62        "with fp = s. Probable cause: s too small. (abs(fp-s)/s>0.001)",
63        ValueError],
64    3: ["The maximal number of iterations (20) allowed for finding smoothing\n"
65        "spline with fp=s has been reached. Probable cause: s too small.\n"
66        "(abs(fp-s)/s>0.001)", ValueError],
67    10: ["Error on input data", ValueError],
68    'unknown': ["An error occurred", TypeError]
69}
70
71_iermess2 = {
72    0: ["The spline has a residual sum of squares fp such that "
73        "abs(fp-s)/s<=0.001", None],
74    -1: ["The spline is an interpolating spline (fp=0)", None],
75    -2: ["The spline is weighted least-squares polynomial of degree kx and ky."
76         "\nfp gives the upper bound fp0 for the smoothing factor s", None],
77    -3: ["Warning. The coefficients of the spline have been computed as the\n"
78         "minimal norm least-squares solution of a rank deficient system.",
79         None],
80    1: ["The required storage space exceeds the available storage space.\n"
81        "Probable causes: nxest or nyest too small or s is too small. (fp>s)",
82        ValueError],
83    2: ["A theoretically impossible result when finding a smoothing spline\n"
84        "with fp = s. Probable causes: s too small or badly chosen eps.\n"
85        "(abs(fp-s)/s>0.001)", ValueError],
86    3: ["The maximal number of iterations (20) allowed for finding smoothing\n"
87        "spline with fp=s has been reached. Probable cause: s too small.\n"
88        "(abs(fp-s)/s>0.001)", ValueError],
89    4: ["No more knots can be added because the number of B-spline\n"
90        "coefficients already exceeds the number of data points m.\n"
91        "Probable causes: either s or m too small. (fp>s)", ValueError],
92    5: ["No more knots can be added because the additional knot would\n"
93        "coincide with an old one. Probable cause: s too small or too large\n"
94        "a weight to an inaccurate data point. (fp>s)", ValueError],
95    10: ["Error on input data", ValueError],
96    11: ["rwrk2 too small, i.e., there is not enough workspace for computing\n"
97         "the minimal least-squares solution of a rank deficient system of\n"
98         "linear equations.", ValueError],
99    'unknown': ["An error occurred", TypeError]
100}
101
102_parcur_cache = {'t': array([], float), 'wrk': array([], float),
103                 'iwrk': array([], dfitpack_int), 'u': array([], float),
104                 'ub': 0, 'ue': 1}
105
106
107def splprep(x, w=None, u=None, ub=None, ue=None, k=3, task=0, s=None, t=None,
108            full_output=0, nest=None, per=0, quiet=1):
109    """
110    Find the B-spline representation of an N-D curve.
111
112    Given a list of N rank-1 arrays, `x`, which represent a curve in
113    N-dimensional space parametrized by `u`, find a smooth approximating
114    spline curve g(`u`). Uses the FORTRAN routine parcur from FITPACK.
115
116    Parameters
117    ----------
118    x : array_like
119        A list of sample vector arrays representing the curve.
120    w : array_like, optional
121        Strictly positive rank-1 array of weights the same length as `x[0]`.
122        The weights are used in computing the weighted least-squares spline
123        fit. If the errors in the `x` values have standard-deviation given by
124        the vector d, then `w` should be 1/d. Default is ``ones(len(x[0]))``.
125    u : array_like, optional
126        An array of parameter values. If not given, these values are
127        calculated automatically as ``M = len(x[0])``, where
128
129            v[0] = 0
130
131            v[i] = v[i-1] + distance(`x[i]`, `x[i-1]`)
132
133            u[i] = v[i] / v[M-1]
134
135    ub, ue : int, optional
136        The end-points of the parameters interval. Defaults to
137        u[0] and u[-1].
138    k : int, optional
139        Degree of the spline. Cubic splines are recommended.
140        Even values of `k` should be avoided especially with a small s-value.
141        ``1 <= k <= 5``, default is 3.
142    task : int, optional
143        If task==0 (default), find t and c for a given smoothing factor, s.
144        If task==1, find t and c for another value of the smoothing factor, s.
145        There must have been a previous call with task=0 or task=1
146        for the same set of data.
147        If task=-1 find the weighted least square spline for a given set of
148        knots, t.
149    s : float, optional
150        A smoothing condition. The amount of smoothness is determined by
151        satisfying the conditions: ``sum((w * (y - g))**2,axis=0) <= s``,
152        where g(x) is the smoothed interpolation of (x,y).  The user can
153        use `s` to control the trade-off between closeness and smoothness
154        of fit. Larger `s` means more smoothing while smaller values of `s`
155        indicate less smoothing. Recommended values of `s` depend on the
156        weights, w.  If the weights represent the inverse of the
157        standard-deviation of y, then a good `s` value should be found in
158        the range ``(m-sqrt(2*m),m+sqrt(2*m))``, where m is the number of
159        data points in x, y, and w.
160    t : int, optional
161        The knots needed for task=-1.
162    full_output : int, optional
163        If non-zero, then return optional outputs.
164    nest : int, optional
165        An over-estimate of the total number of knots of the spline to
166        help in determining the storage space. By default nest=m/2.
167        Always large enough is nest=m+k+1.
168    per : int, optional
169       If non-zero, data points are considered periodic with period
170       ``x[m-1] - x[0]`` and a smooth periodic spline approximation is
171       returned.  Values of ``y[m-1]`` and ``w[m-1]`` are not used.
172    quiet : int, optional
173         Non-zero to suppress messages.
174         This parameter is deprecated; use standard Python warning filters
175         instead.
176
177    Returns
178    -------
179    tck : tuple
180        A tuple (t,c,k) containing the vector of knots, the B-spline
181        coefficients, and the degree of the spline.
182    u : array
183        An array of the values of the parameter.
184    fp : float
185        The weighted sum of squared residuals of the spline approximation.
186    ier : int
187        An integer flag about splrep success.  Success is indicated
188        if ier<=0. If ier in [1,2,3] an error occurred but was not raised.
189        Otherwise an error is raised.
190    msg : str
191        A message corresponding to the integer flag, ier.
192
193    See Also
194    --------
195    splrep, splev, sproot, spalde, splint,
196    bisplrep, bisplev
197    UnivariateSpline, BivariateSpline
198
199    Notes
200    -----
201    See `splev` for evaluation of the spline and its derivatives.
202    The number of dimensions N must be smaller than 11.
203
204    References
205    ----------
206    .. [1] P. Dierckx, "Algorithms for smoothing data with periodic and
207        parametric splines, Computer Graphics and Image Processing",
208        20 (1982) 171-184.
209    .. [2] P. Dierckx, "Algorithms for smoothing data with periodic and
210        parametric splines", report tw55, Dept. Computer Science,
211        K.U.Leuven, 1981.
212    .. [3] P. Dierckx, "Curve and surface fitting with splines", Monographs on
213        Numerical Analysis, Oxford University Press, 1993.
214
215    """
216    if task <= 0:
217        _parcur_cache = {'t': array([], float), 'wrk': array([], float),
218                         'iwrk': array([], dfitpack_int), 'u': array([], float),
219                         'ub': 0, 'ue': 1}
220    x = atleast_1d(x)
221    idim, m = x.shape
222    if per:
223        for i in range(idim):
224            if x[i][0] != x[i][-1]:
225                if quiet < 2:
226                    warnings.warn(RuntimeWarning('Setting x[%d][%d]=x[%d][0]' %
227                                                 (i, m, i)))
228                x[i][-1] = x[i][0]
229    if not 0 < idim < 11:
230        raise TypeError('0 < idim < 11 must hold')
231    if w is None:
232        w = ones(m, float)
233    else:
234        w = atleast_1d(w)
235    ipar = (u is not None)
236    if ipar:
237        _parcur_cache['u'] = u
238        if ub is None:
239            _parcur_cache['ub'] = u[0]
240        else:
241            _parcur_cache['ub'] = ub
242        if ue is None:
243            _parcur_cache['ue'] = u[-1]
244        else:
245            _parcur_cache['ue'] = ue
246    else:
247        _parcur_cache['u'] = zeros(m, float)
248    if not (1 <= k <= 5):
249        raise TypeError('1 <= k= %d <=5 must hold' % k)
250    if not (-1 <= task <= 1):
251        raise TypeError('task must be -1, 0 or 1')
252    if (not len(w) == m) or (ipar == 1 and (not len(u) == m)):
253        raise TypeError('Mismatch of input dimensions')
254    if s is None:
255        s = m - sqrt(2*m)
256    if t is None and task == -1:
257        raise TypeError('Knots must be given for task=-1')
258    if t is not None:
259        _parcur_cache['t'] = atleast_1d(t)
260    n = len(_parcur_cache['t'])
261    if task == -1 and n < 2*k + 2:
262        raise TypeError('There must be at least 2*k+2 knots for task=-1')
263    if m <= k:
264        raise TypeError('m > k must hold')
265    if nest is None:
266        nest = m + 2*k
267
268    if (task >= 0 and s == 0) or (nest < 0):
269        if per:
270            nest = m + 2*k
271        else:
272            nest = m + k + 1
273    nest = max(nest, 2*k + 3)
274    u = _parcur_cache['u']
275    ub = _parcur_cache['ub']
276    ue = _parcur_cache['ue']
277    t = _parcur_cache['t']
278    wrk = _parcur_cache['wrk']
279    iwrk = _parcur_cache['iwrk']
280    t, c, o = _fitpack._parcur(ravel(transpose(x)), w, u, ub, ue, k,
281                               task, ipar, s, t, nest, wrk, iwrk, per)
282    _parcur_cache['u'] = o['u']
283    _parcur_cache['ub'] = o['ub']
284    _parcur_cache['ue'] = o['ue']
285    _parcur_cache['t'] = t
286    _parcur_cache['wrk'] = o['wrk']
287    _parcur_cache['iwrk'] = o['iwrk']
288    ier = o['ier']
289    fp = o['fp']
290    n = len(t)
291    u = o['u']
292    c.shape = idim, n - k - 1
293    tcku = [t, list(c), k], u
294    if ier <= 0 and not quiet:
295        warnings.warn(RuntimeWarning(_iermess[ier][0] +
296                                     "\tk=%d n=%d m=%d fp=%f s=%f" %
297                                     (k, len(t), m, fp, s)))
298    if ier > 0 and not full_output:
299        if ier in [1, 2, 3]:
300            warnings.warn(RuntimeWarning(_iermess[ier][0]))
301        else:
302            try:
303                raise _iermess[ier][1](_iermess[ier][0])
304            except KeyError as e:
305                raise _iermess['unknown'][1](_iermess['unknown'][0]) from e
306    if full_output:
307        try:
308            return tcku, fp, ier, _iermess[ier][0]
309        except KeyError:
310            return tcku, fp, ier, _iermess['unknown'][0]
311    else:
312        return tcku
313
314
315_curfit_cache = {'t': array([], float), 'wrk': array([], float),
316                 'iwrk': array([], dfitpack_int)}
317
318
319def splrep(x, y, w=None, xb=None, xe=None, k=3, task=0, s=None, t=None,
320           full_output=0, per=0, quiet=1):
321    """
322    Find the B-spline representation of 1-D curve.
323
324    Given the set of data points ``(x[i], y[i])`` determine a smooth spline
325    approximation of degree k on the interval ``xb <= x <= xe``.
326
327    Parameters
328    ----------
329    x, y : array_like
330        The data points defining a curve y = f(x).
331    w : array_like, optional
332        Strictly positive rank-1 array of weights the same length as x and y.
333        The weights are used in computing the weighted least-squares spline
334        fit. If the errors in the y values have standard-deviation given by the
335        vector d, then w should be 1/d. Default is ones(len(x)).
336    xb, xe : float, optional
337        The interval to fit.  If None, these default to x[0] and x[-1]
338        respectively.
339    k : int, optional
340        The order of the spline fit. It is recommended to use cubic splines.
341        Even order splines should be avoided especially with small s values.
342        1 <= k <= 5
343    task : {1, 0, -1}, optional
344        If task==0 find t and c for a given smoothing factor, s.
345
346        If task==1 find t and c for another value of the smoothing factor, s.
347        There must have been a previous call with task=0 or task=1 for the same
348        set of data (t will be stored an used internally)
349
350        If task=-1 find the weighted least square spline for a given set of
351        knots, t. These should be interior knots as knots on the ends will be
352        added automatically.
353    s : float, optional
354        A smoothing condition. The amount of smoothness is determined by
355        satisfying the conditions: sum((w * (y - g))**2,axis=0) <= s, where g(x)
356        is the smoothed interpolation of (x,y). The user can use s to control
357        the tradeoff between closeness and smoothness of fit. Larger s means
358        more smoothing while smaller values of s indicate less smoothing.
359        Recommended values of s depend on the weights, w. If the weights
360        represent the inverse of the standard-deviation of y, then a good s
361        value should be found in the range (m-sqrt(2*m),m+sqrt(2*m)) where m is
362        the number of datapoints in x, y, and w. default : s=m-sqrt(2*m) if
363        weights are supplied. s = 0.0 (interpolating) if no weights are
364        supplied.
365    t : array_like, optional
366        The knots needed for task=-1. If given then task is automatically set
367        to -1.
368    full_output : bool, optional
369        If non-zero, then return optional outputs.
370    per : bool, optional
371        If non-zero, data points are considered periodic with period x[m-1] -
372        x[0] and a smooth periodic spline approximation is returned. Values of
373        y[m-1] and w[m-1] are not used.
374    quiet : bool, optional
375        Non-zero to suppress messages.
376        This parameter is deprecated; use standard Python warning filters
377        instead.
378
379    Returns
380    -------
381    tck : tuple
382        (t,c,k) a tuple containing the vector of knots, the B-spline
383        coefficients, and the degree of the spline.
384    fp : array, optional
385        The weighted sum of squared residuals of the spline approximation.
386    ier : int, optional
387        An integer flag about splrep success. Success is indicated if ier<=0.
388        If ier in [1,2,3] an error occurred but was not raised. Otherwise an
389        error is raised.
390    msg : str, optional
391        A message corresponding to the integer flag, ier.
392
393    See Also
394    --------
395    UnivariateSpline, BivariateSpline
396    splprep, splev, sproot, spalde, splint
397    bisplrep, bisplev
398
399    Notes
400    -----
401    See splev for evaluation of the spline and its derivatives. Uses the
402    FORTRAN routine curfit from FITPACK.
403
404    The user is responsible for assuring that the values of *x* are unique.
405    Otherwise, *splrep* will not return sensible results.
406
407    If provided, knots `t` must satisfy the Schoenberg-Whitney conditions,
408    i.e., there must be a subset of data points ``x[j]`` such that
409    ``t[j] < x[j] < t[j+k+1]``, for ``j=0, 1,...,n-k-2``.
410
411    References
412    ----------
413    Based on algorithms described in [1]_, [2]_, [3]_, and [4]_:
414
415    .. [1] P. Dierckx, "An algorithm for smoothing, differentiation and
416       integration of experimental data using spline functions",
417       J.Comp.Appl.Maths 1 (1975) 165-184.
418    .. [2] P. Dierckx, "A fast algorithm for smoothing data on a rectangular
419       grid while using spline functions", SIAM J.Numer.Anal. 19 (1982)
420       1286-1304.
421    .. [3] P. Dierckx, "An improved algorithm for curve fitting with spline
422       functions", report tw54, Dept. Computer Science,K.U. Leuven, 1981.
423    .. [4] P. Dierckx, "Curve and surface fitting with splines", Monographs on
424       Numerical Analysis, Oxford University Press, 1993.
425
426    Examples
427    --------
428
429    >>> import matplotlib.pyplot as plt
430    >>> from scipy.interpolate import splev, splrep
431    >>> x = np.linspace(0, 10, 10)
432    >>> y = np.sin(x)
433    >>> tck = splrep(x, y)
434    >>> x2 = np.linspace(0, 10, 200)
435    >>> y2 = splev(x2, tck)
436    >>> plt.plot(x, y, 'o', x2, y2)
437    >>> plt.show()
438
439    """
440    if task <= 0:
441        _curfit_cache = {}
442    x, y = map(atleast_1d, [x, y])
443    m = len(x)
444    if w is None:
445        w = ones(m, float)
446        if s is None:
447            s = 0.0
448    else:
449        w = atleast_1d(w)
450        if s is None:
451            s = m - sqrt(2*m)
452    if not len(w) == m:
453        raise TypeError('len(w)=%d is not equal to m=%d' % (len(w), m))
454    if (m != len(y)) or (m != len(w)):
455        raise TypeError('Lengths of the first three arguments (x,y,w) must '
456                        'be equal')
457    if not (1 <= k <= 5):
458        raise TypeError('Given degree of the spline (k=%d) is not supported. '
459                        '(1<=k<=5)' % k)
460    if m <= k:
461        raise TypeError('m > k must hold')
462    if xb is None:
463        xb = x[0]
464    if xe is None:
465        xe = x[-1]
466    if not (-1 <= task <= 1):
467        raise TypeError('task must be -1, 0 or 1')
468    if t is not None:
469        task = -1
470    if task == -1:
471        if t is None:
472            raise TypeError('Knots must be given for task=-1')
473        numknots = len(t)
474        _curfit_cache['t'] = empty((numknots + 2*k + 2,), float)
475        _curfit_cache['t'][k+1:-k-1] = t
476        nest = len(_curfit_cache['t'])
477    elif task == 0:
478        if per:
479            nest = max(m + 2*k, 2*k + 3)
480        else:
481            nest = max(m + k + 1, 2*k + 3)
482        t = empty((nest,), float)
483        _curfit_cache['t'] = t
484    if task <= 0:
485        if per:
486            _curfit_cache['wrk'] = empty((m*(k + 1) + nest*(8 + 5*k),), float)
487        else:
488            _curfit_cache['wrk'] = empty((m*(k + 1) + nest*(7 + 3*k),), float)
489        _curfit_cache['iwrk'] = empty((nest,), dfitpack_int)
490    try:
491        t = _curfit_cache['t']
492        wrk = _curfit_cache['wrk']
493        iwrk = _curfit_cache['iwrk']
494    except KeyError as e:
495        raise TypeError("must call with task=1 only after"
496                        " call with task=0,-1") from e
497    if not per:
498        n, c, fp, ier = dfitpack.curfit(task, x, y, w, t, wrk, iwrk,
499                                        xb, xe, k, s)
500    else:
501        n, c, fp, ier = dfitpack.percur(task, x, y, w, t, wrk, iwrk, k, s)
502    tck = (t[:n], c[:n], k)
503    if ier <= 0 and not quiet:
504        _mess = (_iermess[ier][0] + "\tk=%d n=%d m=%d fp=%f s=%f" %
505                 (k, len(t), m, fp, s))
506        warnings.warn(RuntimeWarning(_mess))
507    if ier > 0 and not full_output:
508        if ier in [1, 2, 3]:
509            warnings.warn(RuntimeWarning(_iermess[ier][0]))
510        else:
511            try:
512                raise _iermess[ier][1](_iermess[ier][0])
513            except KeyError as e:
514                raise _iermess['unknown'][1](_iermess['unknown'][0]) from e
515    if full_output:
516        try:
517            return tck, fp, ier, _iermess[ier][0]
518        except KeyError:
519            return tck, fp, ier, _iermess['unknown'][0]
520    else:
521        return tck
522
523
524def splev(x, tck, der=0, ext=0):
525    """
526    Evaluate a B-spline or its derivatives.
527
528    Given the knots and coefficients of a B-spline representation, evaluate
529    the value of the smoothing polynomial and its derivatives. This is a
530    wrapper around the FORTRAN routines splev and splder of FITPACK.
531
532    Parameters
533    ----------
534    x : array_like
535        An array of points at which to return the value of the smoothed
536        spline or its derivatives. If `tck` was returned from `splprep`,
537        then the parameter values, u should be given.
538    tck : tuple
539        A sequence of length 3 returned by `splrep` or `splprep` containing
540        the knots, coefficients, and degree of the spline.
541    der : int, optional
542        The order of derivative of the spline to compute (must be less than
543        or equal to k).
544    ext : int, optional
545        Controls the value returned for elements of ``x`` not in the
546        interval defined by the knot sequence.
547
548        * if ext=0, return the extrapolated value.
549        * if ext=1, return 0
550        * if ext=2, raise a ValueError
551        * if ext=3, return the boundary value.
552
553        The default value is 0.
554
555    Returns
556    -------
557    y : ndarray or list of ndarrays
558        An array of values representing the spline function evaluated at
559        the points in ``x``.  If `tck` was returned from `splprep`, then this
560        is a list of arrays representing the curve in N-D space.
561
562    See Also
563    --------
564    splprep, splrep, sproot, spalde, splint
565    bisplrep, bisplev
566
567    References
568    ----------
569    .. [1] C. de Boor, "On calculating with b-splines", J. Approximation
570        Theory, 6, p.50-62, 1972.
571    .. [2] M.G. Cox, "The numerical evaluation of b-splines", J. Inst. Maths
572        Applics, 10, p.134-149, 1972.
573    .. [3] P. Dierckx, "Curve and surface fitting with splines", Monographs
574        on Numerical Analysis, Oxford University Press, 1993.
575
576    """
577    t, c, k = tck
578    try:
579        c[0][0]
580        parametric = True
581    except Exception:
582        parametric = False
583    if parametric:
584        return list(map(lambda c, x=x, t=t, k=k, der=der:
585                        splev(x, [t, c, k], der, ext), c))
586    else:
587        if not (0 <= der <= k):
588            raise ValueError("0<=der=%d<=k=%d must hold" % (der, k))
589        if ext not in (0, 1, 2, 3):
590            raise ValueError("ext = %s not in (0, 1, 2, 3) " % ext)
591
592        x = asarray(x)
593        shape = x.shape
594        x = atleast_1d(x).ravel()
595        y, ier = _fitpack._spl_(x, der, t, c, k, ext)
596
597        if ier == 10:
598            raise ValueError("Invalid input data")
599        if ier == 1:
600            raise ValueError("Found x value not in the domain")
601        if ier:
602            raise TypeError("An error occurred")
603
604        return y.reshape(shape)
605
606
607def splint(a, b, tck, full_output=0):
608    """
609    Evaluate the definite integral of a B-spline.
610
611    Given the knots and coefficients of a B-spline, evaluate the definite
612    integral of the smoothing polynomial between two given points.
613
614    Parameters
615    ----------
616    a, b : float
617        The end-points of the integration interval.
618    tck : tuple
619        A tuple (t,c,k) containing the vector of knots, the B-spline
620        coefficients, and the degree of the spline (see `splev`).
621    full_output : int, optional
622        Non-zero to return optional output.
623
624    Returns
625    -------
626    integral : float
627        The resulting integral.
628    wrk : ndarray
629        An array containing the integrals of the normalized B-splines
630        defined on the set of knots.
631
632    Notes
633    -----
634    splint silently assumes that the spline function is zero outside the data
635    interval (a, b).
636
637    See Also
638    --------
639    splprep, splrep, sproot, spalde, splev
640    bisplrep, bisplev
641    UnivariateSpline, BivariateSpline
642
643    References
644    ----------
645    .. [1] P.W. Gaffney, The calculation of indefinite integrals of b-splines",
646        J. Inst. Maths Applics, 17, p.37-41, 1976.
647    .. [2] P. Dierckx, "Curve and surface fitting with splines", Monographs
648        on Numerical Analysis, Oxford University Press, 1993.
649
650    """
651    t, c, k = tck
652    try:
653        c[0][0]
654        parametric = True
655    except Exception:
656        parametric = False
657    if parametric:
658        return list(map(lambda c, a=a, b=b, t=t, k=k:
659                        splint(a, b, [t, c, k]), c))
660    else:
661        aint, wrk = _fitpack._splint(t, c, k, a, b)
662        if full_output:
663            return aint, wrk
664        else:
665            return aint
666
667
668def sproot(tck, mest=10):
669    """
670    Find the roots of a cubic B-spline.
671
672    Given the knots (>=8) and coefficients of a cubic B-spline return the
673    roots of the spline.
674
675    Parameters
676    ----------
677    tck : tuple
678        A tuple (t,c,k) containing the vector of knots,
679        the B-spline coefficients, and the degree of the spline.
680        The number of knots must be >= 8, and the degree must be 3.
681        The knots must be a montonically increasing sequence.
682    mest : int, optional
683        An estimate of the number of zeros (Default is 10).
684
685    Returns
686    -------
687    zeros : ndarray
688        An array giving the roots of the spline.
689
690    See also
691    --------
692    splprep, splrep, splint, spalde, splev
693    bisplrep, bisplev
694    UnivariateSpline, BivariateSpline
695
696
697    References
698    ----------
699    .. [1] C. de Boor, "On calculating with b-splines", J. Approximation
700        Theory, 6, p.50-62, 1972.
701    .. [2] M.G. Cox, "The numerical evaluation of b-splines", J. Inst. Maths
702        Applics, 10, p.134-149, 1972.
703    .. [3] P. Dierckx, "Curve and surface fitting with splines", Monographs
704        on Numerical Analysis, Oxford University Press, 1993.
705
706    """
707    t, c, k = tck
708    if k != 3:
709        raise ValueError("sproot works only for cubic (k=3) splines")
710    try:
711        c[0][0]
712        parametric = True
713    except Exception:
714        parametric = False
715    if parametric:
716        return list(map(lambda c, t=t, k=k, mest=mest:
717                        sproot([t, c, k], mest), c))
718    else:
719        if len(t) < 8:
720            raise TypeError("The number of knots %d>=8" % len(t))
721        z, ier = _fitpack._sproot(t, c, k, mest)
722        if ier == 10:
723            raise TypeError("Invalid input data. "
724                            "t1<=..<=t4<t5<..<tn-3<=..<=tn must hold.")
725        if ier == 0:
726            return z
727        if ier == 1:
728            warnings.warn(RuntimeWarning("The number of zeros exceeds mest"))
729            return z
730        raise TypeError("Unknown error")
731
732
733def spalde(x, tck):
734    """
735    Evaluate all derivatives of a B-spline.
736
737    Given the knots and coefficients of a cubic B-spline compute all
738    derivatives up to order k at a point (or set of points).
739
740    Parameters
741    ----------
742    x : array_like
743        A point or a set of points at which to evaluate the derivatives.
744        Note that ``t(k) <= x <= t(n-k+1)`` must hold for each `x`.
745    tck : tuple
746        A tuple (t,c,k) containing the vector of knots,
747        the B-spline coefficients, and the degree of the spline.
748
749    Returns
750    -------
751    results : {ndarray, list of ndarrays}
752        An array (or a list of arrays) containing all derivatives
753        up to order k inclusive for each point `x`.
754
755    See Also
756    --------
757    splprep, splrep, splint, sproot, splev, bisplrep, bisplev,
758    UnivariateSpline, BivariateSpline
759
760    References
761    ----------
762    .. [1] de Boor C : On calculating with b-splines, J. Approximation Theory
763       6 (1972) 50-62.
764    .. [2] Cox M.G. : The numerical evaluation of b-splines, J. Inst. Maths
765       applics 10 (1972) 134-149.
766    .. [3] Dierckx P. : Curve and surface fitting with splines, Monographs on
767       Numerical Analysis, Oxford University Press, 1993.
768
769    """
770    t, c, k = tck
771    try:
772        c[0][0]
773        parametric = True
774    except Exception:
775        parametric = False
776    if parametric:
777        return list(map(lambda c, x=x, t=t, k=k:
778                        spalde(x, [t, c, k]), c))
779    else:
780        x = atleast_1d(x)
781        if len(x) > 1:
782            return list(map(lambda x, tck=tck: spalde(x, tck), x))
783        d, ier = _fitpack._spalde(t, c, k, x[0])
784        if ier == 0:
785            return d
786        if ier == 10:
787            raise TypeError("Invalid input data. t(k)<=x<=t(n-k+1) must hold.")
788        raise TypeError("Unknown error")
789
790# def _curfit(x,y,w=None,xb=None,xe=None,k=3,task=0,s=None,t=None,
791#           full_output=0,nest=None,per=0,quiet=1):
792
793
794_surfit_cache = {'tx': array([], float), 'ty': array([], float),
795                 'wrk': array([], float), 'iwrk': array([], dfitpack_int)}
796
797
798def bisplrep(x, y, z, w=None, xb=None, xe=None, yb=None, ye=None,
799             kx=3, ky=3, task=0, s=None, eps=1e-16, tx=None, ty=None,
800             full_output=0, nxest=None, nyest=None, quiet=1):
801    """
802    Find a bivariate B-spline representation of a surface.
803
804    Given a set of data points (x[i], y[i], z[i]) representing a surface
805    z=f(x,y), compute a B-spline representation of the surface. Based on
806    the routine SURFIT from FITPACK.
807
808    Parameters
809    ----------
810    x, y, z : ndarray
811        Rank-1 arrays of data points.
812    w : ndarray, optional
813        Rank-1 array of weights. By default ``w=np.ones(len(x))``.
814    xb, xe : float, optional
815        End points of approximation interval in `x`.
816        By default ``xb = x.min(), xe=x.max()``.
817    yb, ye : float, optional
818        End points of approximation interval in `y`.
819        By default ``yb=y.min(), ye = y.max()``.
820    kx, ky : int, optional
821        The degrees of the spline (1 <= kx, ky <= 5).
822        Third order (kx=ky=3) is recommended.
823    task : int, optional
824        If task=0, find knots in x and y and coefficients for a given
825        smoothing factor, s.
826        If task=1, find knots and coefficients for another value of the
827        smoothing factor, s.  bisplrep must have been previously called
828        with task=0 or task=1.
829        If task=-1, find coefficients for a given set of knots tx, ty.
830    s : float, optional
831        A non-negative smoothing factor. If weights correspond
832        to the inverse of the standard-deviation of the errors in z,
833        then a good s-value should be found in the range
834        ``(m-sqrt(2*m),m+sqrt(2*m))`` where m=len(x).
835    eps : float, optional
836        A threshold for determining the effective rank of an
837        over-determined linear system of equations (0 < eps < 1).
838        `eps` is not likely to need changing.
839    tx, ty : ndarray, optional
840        Rank-1 arrays of the knots of the spline for task=-1
841    full_output : int, optional
842        Non-zero to return optional outputs.
843    nxest, nyest : int, optional
844        Over-estimates of the total number of knots. If None then
845        ``nxest = max(kx+sqrt(m/2),2*kx+3)``,
846        ``nyest = max(ky+sqrt(m/2),2*ky+3)``.
847    quiet : int, optional
848        Non-zero to suppress printing of messages.
849        This parameter is deprecated; use standard Python warning filters
850        instead.
851
852    Returns
853    -------
854    tck : array_like
855        A list [tx, ty, c, kx, ky] containing the knots (tx, ty) and
856        coefficients (c) of the bivariate B-spline representation of the
857        surface along with the degree of the spline.
858    fp : ndarray
859        The weighted sum of squared residuals of the spline approximation.
860    ier : int
861        An integer flag about splrep success. Success is indicated if
862        ier<=0. If ier in [1,2,3] an error occurred but was not raised.
863        Otherwise an error is raised.
864    msg : str
865        A message corresponding to the integer flag, ier.
866
867    See Also
868    --------
869    splprep, splrep, splint, sproot, splev
870    UnivariateSpline, BivariateSpline
871
872    Notes
873    -----
874    See `bisplev` to evaluate the value of the B-spline given its tck
875    representation.
876
877    References
878    ----------
879    .. [1] Dierckx P.:An algorithm for surface fitting with spline functions
880       Ima J. Numer. Anal. 1 (1981) 267-283.
881    .. [2] Dierckx P.:An algorithm for surface fitting with spline functions
882       report tw50, Dept. Computer Science,K.U.Leuven, 1980.
883    .. [3] Dierckx P.:Curve and surface fitting with splines, Monographs on
884       Numerical Analysis, Oxford University Press, 1993.
885
886    Examples
887    --------
888    Examples are given :ref:`in the tutorial <tutorial-interpolate_2d_spline>`.
889
890    """
891    x, y, z = map(ravel, [x, y, z])  # ensure 1-d arrays.
892    m = len(x)
893    if not (m == len(y) == len(z)):
894        raise TypeError('len(x)==len(y)==len(z) must hold.')
895    if w is None:
896        w = ones(m, float)
897    else:
898        w = atleast_1d(w)
899    if not len(w) == m:
900        raise TypeError('len(w)=%d is not equal to m=%d' % (len(w), m))
901    if xb is None:
902        xb = x.min()
903    if xe is None:
904        xe = x.max()
905    if yb is None:
906        yb = y.min()
907    if ye is None:
908        ye = y.max()
909    if not (-1 <= task <= 1):
910        raise TypeError('task must be -1, 0 or 1')
911    if s is None:
912        s = m - sqrt(2*m)
913    if tx is None and task == -1:
914        raise TypeError('Knots_x must be given for task=-1')
915    if tx is not None:
916        _surfit_cache['tx'] = atleast_1d(tx)
917    nx = len(_surfit_cache['tx'])
918    if ty is None and task == -1:
919        raise TypeError('Knots_y must be given for task=-1')
920    if ty is not None:
921        _surfit_cache['ty'] = atleast_1d(ty)
922    ny = len(_surfit_cache['ty'])
923    if task == -1 and nx < 2*kx+2:
924        raise TypeError('There must be at least 2*kx+2 knots_x for task=-1')
925    if task == -1 and ny < 2*ky+2:
926        raise TypeError('There must be at least 2*ky+2 knots_x for task=-1')
927    if not ((1 <= kx <= 5) and (1 <= ky <= 5)):
928        raise TypeError('Given degree of the spline (kx,ky=%d,%d) is not '
929                        'supported. (1<=k<=5)' % (kx, ky))
930    if m < (kx + 1)*(ky + 1):
931        raise TypeError('m >= (kx+1)(ky+1) must hold')
932    if nxest is None:
933        nxest = int(kx + sqrt(m/2))
934    if nyest is None:
935        nyest = int(ky + sqrt(m/2))
936    nxest, nyest = max(nxest, 2*kx + 3), max(nyest, 2*ky + 3)
937    if task >= 0 and s == 0:
938        nxest = int(kx + sqrt(3*m))
939        nyest = int(ky + sqrt(3*m))
940    if task == -1:
941        _surfit_cache['tx'] = atleast_1d(tx)
942        _surfit_cache['ty'] = atleast_1d(ty)
943    tx, ty = _surfit_cache['tx'], _surfit_cache['ty']
944    wrk = _surfit_cache['wrk']
945    u = nxest - kx - 1
946    v = nyest - ky - 1
947    km = max(kx, ky) + 1
948    ne = max(nxest, nyest)
949    bx, by = kx*v + ky + 1, ky*u + kx + 1
950    b1, b2 = bx, bx + v - ky
951    if bx > by:
952        b1, b2 = by, by + u - kx
953    msg = "Too many data points to interpolate"
954    lwrk1 = _int_overflow(u*v*(2 + b1 + b2) +
955                          2*(u + v + km*(m + ne) + ne - kx - ky) + b2 + 1,
956                          msg=msg)
957    lwrk2 = _int_overflow(u*v*(b2 + 1) + b2, msg=msg)
958    tx, ty, c, o = _fitpack._surfit(x, y, z, w, xb, xe, yb, ye, kx, ky,
959                                    task, s, eps, tx, ty, nxest, nyest,
960                                    wrk, lwrk1, lwrk2)
961    _curfit_cache['tx'] = tx
962    _curfit_cache['ty'] = ty
963    _curfit_cache['wrk'] = o['wrk']
964    ier, fp = o['ier'], o['fp']
965    tck = [tx, ty, c, kx, ky]
966
967    ierm = min(11, max(-3, ier))
968    if ierm <= 0 and not quiet:
969        _mess = (_iermess2[ierm][0] +
970                 "\tkx,ky=%d,%d nx,ny=%d,%d m=%d fp=%f s=%f" %
971                 (kx, ky, len(tx), len(ty), m, fp, s))
972        warnings.warn(RuntimeWarning(_mess))
973    if ierm > 0 and not full_output:
974        if ier in [1, 2, 3, 4, 5]:
975            _mess = ("\n\tkx,ky=%d,%d nx,ny=%d,%d m=%d fp=%f s=%f" %
976                     (kx, ky, len(tx), len(ty), m, fp, s))
977            warnings.warn(RuntimeWarning(_iermess2[ierm][0] + _mess))
978        else:
979            try:
980                raise _iermess2[ierm][1](_iermess2[ierm][0])
981            except KeyError as e:
982                raise _iermess2['unknown'][1](_iermess2['unknown'][0]) from e
983    if full_output:
984        try:
985            return tck, fp, ier, _iermess2[ierm][0]
986        except KeyError:
987            return tck, fp, ier, _iermess2['unknown'][0]
988    else:
989        return tck
990
991
992def bisplev(x, y, tck, dx=0, dy=0):
993    """
994    Evaluate a bivariate B-spline and its derivatives.
995
996    Return a rank-2 array of spline function values (or spline derivative
997    values) at points given by the cross-product of the rank-1 arrays `x` and
998    `y`.  In special cases, return an array or just a float if either `x` or
999    `y` or both are floats.  Based on BISPEV from FITPACK.
1000
1001    Parameters
1002    ----------
1003    x, y : ndarray
1004        Rank-1 arrays specifying the domain over which to evaluate the
1005        spline or its derivative.
1006    tck : tuple
1007        A sequence of length 5 returned by `bisplrep` containing the knot
1008        locations, the coefficients, and the degree of the spline:
1009        [tx, ty, c, kx, ky].
1010    dx, dy : int, optional
1011        The orders of the partial derivatives in `x` and `y` respectively.
1012
1013    Returns
1014    -------
1015    vals : ndarray
1016        The B-spline or its derivative evaluated over the set formed by
1017        the cross-product of `x` and `y`.
1018
1019    See Also
1020    --------
1021    splprep, splrep, splint, sproot, splev
1022    UnivariateSpline, BivariateSpline
1023
1024    Notes
1025    -----
1026        See `bisplrep` to generate the `tck` representation.
1027
1028    References
1029    ----------
1030    .. [1] Dierckx P. : An algorithm for surface fitting
1031       with spline functions
1032       Ima J. Numer. Anal. 1 (1981) 267-283.
1033    .. [2] Dierckx P. : An algorithm for surface fitting
1034       with spline functions
1035       report tw50, Dept. Computer Science,K.U.Leuven, 1980.
1036    .. [3] Dierckx P. : Curve and surface fitting with splines,
1037       Monographs on Numerical Analysis, Oxford University Press, 1993.
1038
1039    Examples
1040    --------
1041    Examples are given :ref:`in the tutorial <tutorial-interpolate_2d_spline>`.
1042
1043    """
1044    tx, ty, c, kx, ky = tck
1045    if not (0 <= dx < kx):
1046        raise ValueError("0 <= dx = %d < kx = %d must hold" % (dx, kx))
1047    if not (0 <= dy < ky):
1048        raise ValueError("0 <= dy = %d < ky = %d must hold" % (dy, ky))
1049    x, y = map(atleast_1d, [x, y])
1050    if (len(x.shape) != 1) or (len(y.shape) != 1):
1051        raise ValueError("First two entries should be rank-1 arrays.")
1052    z, ier = _fitpack._bispev(tx, ty, c, kx, ky, x, y, dx, dy)
1053    if ier == 10:
1054        raise ValueError("Invalid input data")
1055    if ier:
1056        raise TypeError("An error occurred")
1057    z.shape = len(x), len(y)
1058    if len(z) > 1:
1059        return z
1060    if len(z[0]) > 1:
1061        return z[0]
1062    return z[0][0]
1063
1064
1065def dblint(xa, xb, ya, yb, tck):
1066    """Evaluate the integral of a spline over area [xa,xb] x [ya,yb].
1067
1068    Parameters
1069    ----------
1070    xa, xb : float
1071        The end-points of the x integration interval.
1072    ya, yb : float
1073        The end-points of the y integration interval.
1074    tck : list [tx, ty, c, kx, ky]
1075        A sequence of length 5 returned by bisplrep containing the knot
1076        locations tx, ty, the coefficients c, and the degrees kx, ky
1077        of the spline.
1078
1079    Returns
1080    -------
1081    integ : float
1082        The value of the resulting integral.
1083    """
1084    tx, ty, c, kx, ky = tck
1085    return dfitpack.dblint(tx, ty, c, kx, ky, xa, xb, ya, yb)
1086
1087
1088def insert(x, tck, m=1, per=0):
1089    """
1090    Insert knots into a B-spline.
1091
1092    Given the knots and coefficients of a B-spline representation, create a
1093    new B-spline with a knot inserted `m` times at point `x`.
1094    This is a wrapper around the FORTRAN routine insert of FITPACK.
1095
1096    Parameters
1097    ----------
1098    x (u) : array_like
1099        A 1-D point at which to insert a new knot(s).  If `tck` was returned
1100        from ``splprep``, then the parameter values, u should be given.
1101    tck : tuple
1102        A tuple (t,c,k) returned by ``splrep`` or ``splprep`` containing
1103        the vector of knots, the B-spline coefficients,
1104        and the degree of the spline.
1105    m : int, optional
1106        The number of times to insert the given knot (its multiplicity).
1107        Default is 1.
1108    per : int, optional
1109        If non-zero, the input spline is considered periodic.
1110
1111    Returns
1112    -------
1113    tck : tuple
1114        A tuple (t,c,k) containing the vector of knots, the B-spline
1115        coefficients, and the degree of the new spline.
1116        ``t(k+1) <= x <= t(n-k)``, where k is the degree of the spline.
1117        In case of a periodic spline (``per != 0``) there must be
1118        either at least k interior knots t(j) satisfying ``t(k+1)<t(j)<=x``
1119        or at least k interior knots t(j) satisfying ``x<=t(j)<t(n-k)``.
1120
1121    Notes
1122    -----
1123    Based on algorithms from [1]_ and [2]_.
1124
1125    References
1126    ----------
1127    .. [1] W. Boehm, "Inserting new knots into b-spline curves.",
1128        Computer Aided Design, 12, p.199-201, 1980.
1129    .. [2] P. Dierckx, "Curve and surface fitting with splines, Monographs on
1130        Numerical Analysis", Oxford University Press, 1993.
1131
1132    """
1133    t, c, k = tck
1134    try:
1135        c[0][0]
1136        parametric = True
1137    except Exception:
1138        parametric = False
1139    if parametric:
1140        cc = []
1141        for c_vals in c:
1142            tt, cc_val, kk = insert(x, [t, c_vals, k], m)
1143            cc.append(cc_val)
1144        return (tt, cc, kk)
1145    else:
1146        tt, cc, ier = _fitpack._insert(per, t, c, k, x, m)
1147        if ier == 10:
1148            raise ValueError("Invalid input data")
1149        if ier:
1150            raise TypeError("An error occurred")
1151        return (tt, cc, k)
1152
1153
1154def splder(tck, n=1):
1155    """
1156    Compute the spline representation of the derivative of a given spline
1157
1158    Parameters
1159    ----------
1160    tck : tuple of (t, c, k)
1161        Spline whose derivative to compute
1162    n : int, optional
1163        Order of derivative to evaluate. Default: 1
1164
1165    Returns
1166    -------
1167    tck_der : tuple of (t2, c2, k2)
1168        Spline of order k2=k-n representing the derivative
1169        of the input spline.
1170
1171    Notes
1172    -----
1173
1174    .. versionadded:: 0.13.0
1175
1176    See Also
1177    --------
1178    splantider, splev, spalde
1179
1180    Examples
1181    --------
1182    This can be used for finding maxima of a curve:
1183
1184    >>> from scipy.interpolate import splrep, splder, sproot
1185    >>> x = np.linspace(0, 10, 70)
1186    >>> y = np.sin(x)
1187    >>> spl = splrep(x, y, k=4)
1188
1189    Now, differentiate the spline and find the zeros of the
1190    derivative. (NB: `sproot` only works for order 3 splines, so we
1191    fit an order 4 spline):
1192
1193    >>> dspl = splder(spl)
1194    >>> sproot(dspl) / np.pi
1195    array([ 0.50000001,  1.5       ,  2.49999998])
1196
1197    This agrees well with roots :math:`\\pi/2 + n\\pi` of
1198    :math:`\\cos(x) = \\sin'(x)`.
1199
1200    """
1201    if n < 0:
1202        return splantider(tck, -n)
1203
1204    t, c, k = tck
1205
1206    if n > k:
1207        raise ValueError(("Order of derivative (n = %r) must be <= "
1208                          "order of spline (k = %r)") % (n, tck[2]))
1209
1210    # Extra axes for the trailing dims of the `c` array:
1211    sh = (slice(None),) + ((None,)*len(c.shape[1:]))
1212
1213    with np.errstate(invalid='raise', divide='raise'):
1214        try:
1215            for j in range(n):
1216                # See e.g. Schumaker, Spline Functions: Basic Theory, Chapter 5
1217
1218                # Compute the denominator in the differentiation formula.
1219                # (and append traling dims, if necessary)
1220                dt = t[k+1:-1] - t[1:-k-1]
1221                dt = dt[sh]
1222                # Compute the new coefficients
1223                c = (c[1:-1-k] - c[:-2-k]) * k / dt
1224                # Pad coefficient array to same size as knots (FITPACK
1225                # convention)
1226                c = np.r_[c, np.zeros((k,) + c.shape[1:])]
1227                # Adjust knots
1228                t = t[1:-1]
1229                k -= 1
1230        except FloatingPointError as e:
1231            raise ValueError(("The spline has internal repeated knots "
1232                              "and is not differentiable %d times") % n) from e
1233
1234    return t, c, k
1235
1236
1237def splantider(tck, n=1):
1238    """
1239    Compute the spline for the antiderivative (integral) of a given spline.
1240
1241    Parameters
1242    ----------
1243    tck : tuple of (t, c, k)
1244        Spline whose antiderivative to compute
1245    n : int, optional
1246        Order of antiderivative to evaluate. Default: 1
1247
1248    Returns
1249    -------
1250    tck_ader : tuple of (t2, c2, k2)
1251        Spline of order k2=k+n representing the antiderivative of the input
1252        spline.
1253
1254    See Also
1255    --------
1256    splder, splev, spalde
1257
1258    Notes
1259    -----
1260    The `splder` function is the inverse operation of this function.
1261    Namely, ``splder(splantider(tck))`` is identical to `tck`, modulo
1262    rounding error.
1263
1264    .. versionadded:: 0.13.0
1265
1266    Examples
1267    --------
1268    >>> from scipy.interpolate import splrep, splder, splantider, splev
1269    >>> x = np.linspace(0, np.pi/2, 70)
1270    >>> y = 1 / np.sqrt(1 - 0.8*np.sin(x)**2)
1271    >>> spl = splrep(x, y)
1272
1273    The derivative is the inverse operation of the antiderivative,
1274    although some floating point error accumulates:
1275
1276    >>> splev(1.7, spl), splev(1.7, splder(splantider(spl)))
1277    (array(2.1565429877197317), array(2.1565429877201865))
1278
1279    Antiderivative can be used to evaluate definite integrals:
1280
1281    >>> ispl = splantider(spl)
1282    >>> splev(np.pi/2, ispl) - splev(0, ispl)
1283    2.2572053588768486
1284
1285    This is indeed an approximation to the complete elliptic integral
1286    :math:`K(m) = \\int_0^{\\pi/2} [1 - m\\sin^2 x]^{-1/2} dx`:
1287
1288    >>> from scipy.special import ellipk
1289    >>> ellipk(0.8)
1290    2.2572053268208538
1291
1292    """
1293    if n < 0:
1294        return splder(tck, -n)
1295
1296    t, c, k = tck
1297
1298    # Extra axes for the trailing dims of the `c` array:
1299    sh = (slice(None),) + (None,)*len(c.shape[1:])
1300
1301    for j in range(n):
1302        # This is the inverse set of operations to splder.
1303
1304        # Compute the multiplier in the antiderivative formula.
1305        dt = t[k+1:] - t[:-k-1]
1306        dt = dt[sh]
1307        # Compute the new coefficients
1308        c = np.cumsum(c[:-k-1] * dt, axis=0) / (k + 1)
1309        c = np.r_[np.zeros((1,) + c.shape[1:]),
1310                  c,
1311                  [c[-1]] * (k+2)]
1312        # New knots
1313        t = np.r_[t[0], t, t[-1]]
1314        k += 1
1315
1316    return t, c, k
1317