Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

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: 

305 raise _iermess['unknown'][1](_iermess['unknown'][0]) 

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 Notes 

394 ----- 

395 See splev for evaluation of the spline and its derivatives. 

396 

397 The user is responsible for assuring that the values of *x* are unique. 

398 Otherwise, *splrep* will not return sensible results. 

399 

400 See Also 

401 -------- 

402 UnivariateSpline, BivariateSpline 

403 splprep, splev, sproot, spalde, splint 

404 bisplrep, bisplev 

405 

406 Notes 

407 ----- 

408 See splev for evaluation of the spline and its derivatives. Uses the 

409 FORTRAN routine curfit from FITPACK. 

410 

411 If provided, knots `t` must satisfy the Schoenberg-Whitney conditions, 

412 i.e., there must be a subset of data points ``x[j]`` such that 

413 ``t[j] < x[j] < t[j+k+1]``, for ``j=0, 1,...,n-k-2``. 

414 

415 References 

416 ---------- 

417 Based on algorithms described in [1]_, [2]_, [3]_, and [4]_: 

418 

419 .. [1] P. Dierckx, "An algorithm for smoothing, differentiation and 

420 integration of experimental data using spline functions", 

421 J.Comp.Appl.Maths 1 (1975) 165-184. 

422 .. [2] P. Dierckx, "A fast algorithm for smoothing data on a rectangular 

423 grid while using spline functions", SIAM J.Numer.Anal. 19 (1982) 

424 1286-1304. 

425 .. [3] P. Dierckx, "An improved algorithm for curve fitting with spline 

426 functions", report tw54, Dept. Computer Science,K.U. Leuven, 1981. 

427 .. [4] P. Dierckx, "Curve and surface fitting with splines", Monographs on 

428 Numerical Analysis, Oxford University Press, 1993. 

429 

430 Examples 

431 -------- 

432 

433 >>> import matplotlib.pyplot as plt 

434 >>> from scipy.interpolate import splev, splrep 

435 >>> x = np.linspace(0, 10, 10) 

436 >>> y = np.sin(x) 

437 >>> tck = splrep(x, y) 

438 >>> x2 = np.linspace(0, 10, 200) 

439 >>> y2 = splev(x2, tck) 

440 >>> plt.plot(x, y, 'o', x2, y2) 

441 >>> plt.show() 

442 

443 """ 

444 if task <= 0: 

445 _curfit_cache = {} 

446 x, y = map(atleast_1d, [x, y]) 

447 m = len(x) 

448 if w is None: 

449 w = ones(m, float) 

450 if s is None: 

451 s = 0.0 

452 else: 

453 w = atleast_1d(w) 

454 if s is None: 

455 s = m - sqrt(2*m) 

456 if not len(w) == m: 

457 raise TypeError('len(w)=%d is not equal to m=%d' % (len(w), m)) 

458 if (m != len(y)) or (m != len(w)): 

459 raise TypeError('Lengths of the first three arguments (x,y,w) must ' 

460 'be equal') 

461 if not (1 <= k <= 5): 

462 raise TypeError('Given degree of the spline (k=%d) is not supported. ' 

463 '(1<=k<=5)' % k) 

464 if m <= k: 

465 raise TypeError('m > k must hold') 

466 if xb is None: 

467 xb = x[0] 

468 if xe is None: 

469 xe = x[-1] 

470 if not (-1 <= task <= 1): 

471 raise TypeError('task must be -1, 0 or 1') 

472 if t is not None: 

473 task = -1 

474 if task == -1: 

475 if t is None: 

476 raise TypeError('Knots must be given for task=-1') 

477 numknots = len(t) 

478 _curfit_cache['t'] = empty((numknots + 2*k + 2,), float) 

479 _curfit_cache['t'][k+1:-k-1] = t 

480 nest = len(_curfit_cache['t']) 

481 elif task == 0: 

482 if per: 

483 nest = max(m + 2*k, 2*k + 3) 

484 else: 

485 nest = max(m + k + 1, 2*k + 3) 

486 t = empty((nest,), float) 

487 _curfit_cache['t'] = t 

488 if task <= 0: 

489 if per: 

490 _curfit_cache['wrk'] = empty((m*(k + 1) + nest*(8 + 5*k),), float) 

491 else: 

492 _curfit_cache['wrk'] = empty((m*(k + 1) + nest*(7 + 3*k),), float) 

493 _curfit_cache['iwrk'] = empty((nest,), dfitpack_int) 

494 try: 

495 t = _curfit_cache['t'] 

496 wrk = _curfit_cache['wrk'] 

497 iwrk = _curfit_cache['iwrk'] 

498 except KeyError: 

499 raise TypeError("must call with task=1 only after" 

500 " call with task=0,-1") 

501 if not per: 

502 n, c, fp, ier = dfitpack.curfit(task, x, y, w, t, wrk, iwrk, 

503 xb, xe, k, s) 

504 else: 

505 n, c, fp, ier = dfitpack.percur(task, x, y, w, t, wrk, iwrk, k, s) 

506 tck = (t[:n], c[:n], k) 

507 if ier <= 0 and not quiet: 

508 _mess = (_iermess[ier][0] + "\tk=%d n=%d m=%d fp=%f s=%f" % 

509 (k, len(t), m, fp, s)) 

510 warnings.warn(RuntimeWarning(_mess)) 

511 if ier > 0 and not full_output: 

512 if ier in [1, 2, 3]: 

513 warnings.warn(RuntimeWarning(_iermess[ier][0])) 

514 else: 

515 try: 

516 raise _iermess[ier][1](_iermess[ier][0]) 

517 except KeyError: 

518 raise _iermess['unknown'][1](_iermess['unknown'][0]) 

519 if full_output: 

520 try: 

521 return tck, fp, ier, _iermess[ier][0] 

522 except KeyError: 

523 return tck, fp, ier, _iermess['unknown'][0] 

524 else: 

525 return tck 

526 

527 

528def splev(x, tck, der=0, ext=0): 

529 """ 

530 Evaluate a B-spline or its derivatives. 

531 

532 Given the knots and coefficients of a B-spline representation, evaluate 

533 the value of the smoothing polynomial and its derivatives. This is a 

534 wrapper around the FORTRAN routines splev and splder of FITPACK. 

535 

536 Parameters 

537 ---------- 

538 x : array_like 

539 An array of points at which to return the value of the smoothed 

540 spline or its derivatives. If `tck` was returned from `splprep`, 

541 then the parameter values, u should be given. 

542 tck : tuple 

543 A sequence of length 3 returned by `splrep` or `splprep` containing 

544 the knots, coefficients, and degree of the spline. 

545 der : int, optional 

546 The order of derivative of the spline to compute (must be less than 

547 or equal to k). 

548 ext : int, optional 

549 Controls the value returned for elements of ``x`` not in the 

550 interval defined by the knot sequence. 

551 

552 * if ext=0, return the extrapolated value. 

553 * if ext=1, return 0 

554 * if ext=2, raise a ValueError 

555 * if ext=3, return the boundary value. 

556 

557 The default value is 0. 

558 

559 Returns 

560 ------- 

561 y : ndarray or list of ndarrays 

562 An array of values representing the spline function evaluated at 

563 the points in ``x``. If `tck` was returned from `splprep`, then this 

564 is a list of arrays representing the curve in N-D space. 

565 

566 See Also 

567 -------- 

568 splprep, splrep, sproot, spalde, splint 

569 bisplrep, bisplev 

570 

571 References 

572 ---------- 

573 .. [1] C. de Boor, "On calculating with b-splines", J. Approximation 

574 Theory, 6, p.50-62, 1972. 

575 .. [2] M.G. Cox, "The numerical evaluation of b-splines", J. Inst. Maths 

576 Applics, 10, p.134-149, 1972. 

577 .. [3] P. Dierckx, "Curve and surface fitting with splines", Monographs 

578 on Numerical Analysis, Oxford University Press, 1993. 

579 

580 """ 

581 t, c, k = tck 

582 try: 

583 c[0][0] 

584 parametric = True 

585 except Exception: 

586 parametric = False 

587 if parametric: 

588 return list(map(lambda c, x=x, t=t, k=k, der=der: 

589 splev(x, [t, c, k], der, ext), c)) 

590 else: 

591 if not (0 <= der <= k): 

592 raise ValueError("0<=der=%d<=k=%d must hold" % (der, k)) 

593 if ext not in (0, 1, 2, 3): 

594 raise ValueError("ext = %s not in (0, 1, 2, 3) " % ext) 

595 

596 x = asarray(x) 

597 shape = x.shape 

598 x = atleast_1d(x).ravel() 

599 y, ier = _fitpack._spl_(x, der, t, c, k, ext) 

600 

601 if ier == 10: 

602 raise ValueError("Invalid input data") 

603 if ier == 1: 

604 raise ValueError("Found x value not in the domain") 

605 if ier: 

606 raise TypeError("An error occurred") 

607 

608 return y.reshape(shape) 

609 

610 

611def splint(a, b, tck, full_output=0): 

612 """ 

613 Evaluate the definite integral of a B-spline. 

614 

615 Given the knots and coefficients of a B-spline, evaluate the definite 

616 integral of the smoothing polynomial between two given points. 

617 

618 Parameters 

619 ---------- 

620 a, b : float 

621 The end-points of the integration interval. 

622 tck : tuple 

623 A tuple (t,c,k) containing the vector of knots, the B-spline 

624 coefficients, and the degree of the spline (see `splev`). 

625 full_output : int, optional 

626 Non-zero to return optional output. 

627 

628 Returns 

629 ------- 

630 integral : float 

631 The resulting integral. 

632 wrk : ndarray 

633 An array containing the integrals of the normalized B-splines 

634 defined on the set of knots. 

635 

636 Notes 

637 ----- 

638 splint silently assumes that the spline function is zero outside the data 

639 interval (a, b). 

640 

641 See Also 

642 -------- 

643 splprep, splrep, sproot, spalde, splev 

644 bisplrep, bisplev 

645 UnivariateSpline, BivariateSpline 

646 

647 References 

648 ---------- 

649 .. [1] P.W. Gaffney, The calculation of indefinite integrals of b-splines", 

650 J. Inst. Maths Applics, 17, p.37-41, 1976. 

651 .. [2] P. Dierckx, "Curve and surface fitting with splines", Monographs 

652 on Numerical Analysis, Oxford University Press, 1993. 

653 

654 """ 

655 t, c, k = tck 

656 try: 

657 c[0][0] 

658 parametric = True 

659 except Exception: 

660 parametric = False 

661 if parametric: 

662 return list(map(lambda c, a=a, b=b, t=t, k=k: 

663 splint(a, b, [t, c, k]), c)) 

664 else: 

665 aint, wrk = _fitpack._splint(t, c, k, a, b) 

666 if full_output: 

667 return aint, wrk 

668 else: 

669 return aint 

670 

671 

672def sproot(tck, mest=10): 

673 """ 

674 Find the roots of a cubic B-spline. 

675 

676 Given the knots (>=8) and coefficients of a cubic B-spline return the 

677 roots of the spline. 

678 

679 Parameters 

680 ---------- 

681 tck : tuple 

682 A tuple (t,c,k) containing the vector of knots, 

683 the B-spline coefficients, and the degree of the spline. 

684 The number of knots must be >= 8, and the degree must be 3. 

685 The knots must be a montonically increasing sequence. 

686 mest : int, optional 

687 An estimate of the number of zeros (Default is 10). 

688 

689 Returns 

690 ------- 

691 zeros : ndarray 

692 An array giving the roots of the spline. 

693 

694 See also 

695 -------- 

696 splprep, splrep, splint, spalde, splev 

697 bisplrep, bisplev 

698 UnivariateSpline, BivariateSpline 

699 

700 

701 References 

702 ---------- 

703 .. [1] C. de Boor, "On calculating with b-splines", J. Approximation 

704 Theory, 6, p.50-62, 1972. 

705 .. [2] M.G. Cox, "The numerical evaluation of b-splines", J. Inst. Maths 

706 Applics, 10, p.134-149, 1972. 

707 .. [3] P. Dierckx, "Curve and surface fitting with splines", Monographs 

708 on Numerical Analysis, Oxford University Press, 1993. 

709 

710 """ 

711 t, c, k = tck 

712 if k != 3: 

713 raise ValueError("sproot works only for cubic (k=3) splines") 

714 try: 

715 c[0][0] 

716 parametric = True 

717 except Exception: 

718 parametric = False 

719 if parametric: 

720 return list(map(lambda c, t=t, k=k, mest=mest: 

721 sproot([t, c, k], mest), c)) 

722 else: 

723 if len(t) < 8: 

724 raise TypeError("The number of knots %d>=8" % len(t)) 

725 z, ier = _fitpack._sproot(t, c, k, mest) 

726 if ier == 10: 

727 raise TypeError("Invalid input data. " 

728 "t1<=..<=t4<t5<..<tn-3<=..<=tn must hold.") 

729 if ier == 0: 

730 return z 

731 if ier == 1: 

732 warnings.warn(RuntimeWarning("The number of zeros exceeds mest")) 

733 return z 

734 raise TypeError("Unknown error") 

735 

736 

737def spalde(x, tck): 

738 """ 

739 Evaluate all derivatives of a B-spline. 

740 

741 Given the knots and coefficients of a cubic B-spline compute all 

742 derivatives up to order k at a point (or set of points). 

743 

744 Parameters 

745 ---------- 

746 x : array_like 

747 A point or a set of points at which to evaluate the derivatives. 

748 Note that ``t(k) <= x <= t(n-k+1)`` must hold for each `x`. 

749 tck : tuple 

750 A tuple (t,c,k) containing the vector of knots, 

751 the B-spline coefficients, and the degree of the spline. 

752 

753 Returns 

754 ------- 

755 results : {ndarray, list of ndarrays} 

756 An array (or a list of arrays) containing all derivatives 

757 up to order k inclusive for each point `x`. 

758 

759 See Also 

760 -------- 

761 splprep, splrep, splint, sproot, splev, bisplrep, bisplev, 

762 UnivariateSpline, BivariateSpline 

763 

764 References 

765 ---------- 

766 .. [1] de Boor C : On calculating with b-splines, J. Approximation Theory 

767 6 (1972) 50-62. 

768 .. [2] Cox M.G. : The numerical evaluation of b-splines, J. Inst. Maths 

769 applics 10 (1972) 134-149. 

770 .. [3] Dierckx P. : Curve and surface fitting with splines, Monographs on 

771 Numerical Analysis, Oxford University Press, 1993. 

772 

773 """ 

774 t, c, k = tck 

775 try: 

776 c[0][0] 

777 parametric = True 

778 except Exception: 

779 parametric = False 

780 if parametric: 

781 return list(map(lambda c, x=x, t=t, k=k: 

782 spalde(x, [t, c, k]), c)) 

783 else: 

784 x = atleast_1d(x) 

785 if len(x) > 1: 

786 return list(map(lambda x, tck=tck: spalde(x, tck), x)) 

787 d, ier = _fitpack._spalde(t, c, k, x[0]) 

788 if ier == 0: 

789 return d 

790 if ier == 10: 

791 raise TypeError("Invalid input data. t(k)<=x<=t(n-k+1) must hold.") 

792 raise TypeError("Unknown error") 

793 

794# def _curfit(x,y,w=None,xb=None,xe=None,k=3,task=0,s=None,t=None, 

795# full_output=0,nest=None,per=0,quiet=1): 

796 

797 

798_surfit_cache = {'tx': array([], float), 'ty': array([], float), 

799 'wrk': array([], float), 'iwrk': array([], dfitpack_int)} 

800 

801 

802def bisplrep(x, y, z, w=None, xb=None, xe=None, yb=None, ye=None, 

803 kx=3, ky=3, task=0, s=None, eps=1e-16, tx=None, ty=None, 

804 full_output=0, nxest=None, nyest=None, quiet=1): 

805 """ 

806 Find a bivariate B-spline representation of a surface. 

807 

808 Given a set of data points (x[i], y[i], z[i]) representing a surface 

809 z=f(x,y), compute a B-spline representation of the surface. Based on 

810 the routine SURFIT from FITPACK. 

811 

812 Parameters 

813 ---------- 

814 x, y, z : ndarray 

815 Rank-1 arrays of data points. 

816 w : ndarray, optional 

817 Rank-1 array of weights. By default ``w=np.ones(len(x))``. 

818 xb, xe : float, optional 

819 End points of approximation interval in `x`. 

820 By default ``xb = x.min(), xe=x.max()``. 

821 yb, ye : float, optional 

822 End points of approximation interval in `y`. 

823 By default ``yb=y.min(), ye = y.max()``. 

824 kx, ky : int, optional 

825 The degrees of the spline (1 <= kx, ky <= 5). 

826 Third order (kx=ky=3) is recommended. 

827 task : int, optional 

828 If task=0, find knots in x and y and coefficients for a given 

829 smoothing factor, s. 

830 If task=1, find knots and coefficients for another value of the 

831 smoothing factor, s. bisplrep must have been previously called 

832 with task=0 or task=1. 

833 If task=-1, find coefficients for a given set of knots tx, ty. 

834 s : float, optional 

835 A non-negative smoothing factor. If weights correspond 

836 to the inverse of the standard-deviation of the errors in z, 

837 then a good s-value should be found in the range 

838 ``(m-sqrt(2*m),m+sqrt(2*m))`` where m=len(x). 

839 eps : float, optional 

840 A threshold for determining the effective rank of an 

841 over-determined linear system of equations (0 < eps < 1). 

842 `eps` is not likely to need changing. 

843 tx, ty : ndarray, optional 

844 Rank-1 arrays of the knots of the spline for task=-1 

845 full_output : int, optional 

846 Non-zero to return optional outputs. 

847 nxest, nyest : int, optional 

848 Over-estimates of the total number of knots. If None then 

849 ``nxest = max(kx+sqrt(m/2),2*kx+3)``, 

850 ``nyest = max(ky+sqrt(m/2),2*ky+3)``. 

851 quiet : int, optional 

852 Non-zero to suppress printing of messages. 

853 This parameter is deprecated; use standard Python warning filters 

854 instead. 

855 

856 Returns 

857 ------- 

858 tck : array_like 

859 A list [tx, ty, c, kx, ky] containing the knots (tx, ty) and 

860 coefficients (c) of the bivariate B-spline representation of the 

861 surface along with the degree of the spline. 

862 fp : ndarray 

863 The weighted sum of squared residuals of the spline approximation. 

864 ier : int 

865 An integer flag about splrep success. Success is indicated if 

866 ier<=0. If ier in [1,2,3] an error occurred but was not raised. 

867 Otherwise an error is raised. 

868 msg : str 

869 A message corresponding to the integer flag, ier. 

870 

871 See Also 

872 -------- 

873 splprep, splrep, splint, sproot, splev 

874 UnivariateSpline, BivariateSpline 

875 

876 Notes 

877 ----- 

878 See `bisplev` to evaluate the value of the B-spline given its tck 

879 representation. 

880 

881 References 

882 ---------- 

883 .. [1] Dierckx P.:An algorithm for surface fitting with spline functions 

884 Ima J. Numer. Anal. 1 (1981) 267-283. 

885 .. [2] Dierckx P.:An algorithm for surface fitting with spline functions 

886 report tw50, Dept. Computer Science,K.U.Leuven, 1980. 

887 .. [3] Dierckx P.:Curve and surface fitting with splines, Monographs on 

888 Numerical Analysis, Oxford University Press, 1993. 

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: 

982 raise _iermess2['unknown'][1](_iermess2['unknown'][0]) 

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 """ 

1040 tx, ty, c, kx, ky = tck 

1041 if not (0 <= dx < kx): 

1042 raise ValueError("0 <= dx = %d < kx = %d must hold" % (dx, kx)) 

1043 if not (0 <= dy < ky): 

1044 raise ValueError("0 <= dy = %d < ky = %d must hold" % (dy, ky)) 

1045 x, y = map(atleast_1d, [x, y]) 

1046 if (len(x.shape) != 1) or (len(y.shape) != 1): 

1047 raise ValueError("First two entries should be rank-1 arrays.") 

1048 z, ier = _fitpack._bispev(tx, ty, c, kx, ky, x, y, dx, dy) 

1049 if ier == 10: 

1050 raise ValueError("Invalid input data") 

1051 if ier: 

1052 raise TypeError("An error occurred") 

1053 z.shape = len(x), len(y) 

1054 if len(z) > 1: 

1055 return z 

1056 if len(z[0]) > 1: 

1057 return z[0] 

1058 return z[0][0] 

1059 

1060 

1061def dblint(xa, xb, ya, yb, tck): 

1062 """Evaluate the integral of a spline over area [xa,xb] x [ya,yb]. 

1063 

1064 Parameters 

1065 ---------- 

1066 xa, xb : float 

1067 The end-points of the x integration interval. 

1068 ya, yb : float 

1069 The end-points of the y integration interval. 

1070 tck : list [tx, ty, c, kx, ky] 

1071 A sequence of length 5 returned by bisplrep containing the knot 

1072 locations tx, ty, the coefficients c, and the degrees kx, ky 

1073 of the spline. 

1074 

1075 Returns 

1076 ------- 

1077 integ : float 

1078 The value of the resulting integral. 

1079 """ 

1080 tx, ty, c, kx, ky = tck 

1081 return dfitpack.dblint(tx, ty, c, kx, ky, xa, xb, ya, yb) 

1082 

1083 

1084def insert(x, tck, m=1, per=0): 

1085 """ 

1086 Insert knots into a B-spline. 

1087 

1088 Given the knots and coefficients of a B-spline representation, create a 

1089 new B-spline with a knot inserted `m` times at point `x`. 

1090 This is a wrapper around the FORTRAN routine insert of FITPACK. 

1091 

1092 Parameters 

1093 ---------- 

1094 x (u) : array_like 

1095 A 1-D point at which to insert a new knot(s). If `tck` was returned 

1096 from ``splprep``, then the parameter values, u should be given. 

1097 tck : tuple 

1098 A tuple (t,c,k) returned by ``splrep`` or ``splprep`` containing 

1099 the vector of knots, the B-spline coefficients, 

1100 and the degree of the spline. 

1101 m : int, optional 

1102 The number of times to insert the given knot (its multiplicity). 

1103 Default is 1. 

1104 per : int, optional 

1105 If non-zero, the input spline is considered periodic. 

1106 

1107 Returns 

1108 ------- 

1109 tck : tuple 

1110 A tuple (t,c,k) containing the vector of knots, the B-spline 

1111 coefficients, and the degree of the new spline. 

1112 ``t(k+1) <= x <= t(n-k)``, where k is the degree of the spline. 

1113 In case of a periodic spline (``per != 0``) there must be 

1114 either at least k interior knots t(j) satisfying ``t(k+1)<t(j)<=x`` 

1115 or at least k interior knots t(j) satisfying ``x<=t(j)<t(n-k)``. 

1116 

1117 Notes 

1118 ----- 

1119 Based on algorithms from [1]_ and [2]_. 

1120 

1121 References 

1122 ---------- 

1123 .. [1] W. Boehm, "Inserting new knots into b-spline curves.", 

1124 Computer Aided Design, 12, p.199-201, 1980. 

1125 .. [2] P. Dierckx, "Curve and surface fitting with splines, Monographs on 

1126 Numerical Analysis", Oxford University Press, 1993. 

1127 

1128 """ 

1129 t, c, k = tck 

1130 try: 

1131 c[0][0] 

1132 parametric = True 

1133 except Exception: 

1134 parametric = False 

1135 if parametric: 

1136 cc = [] 

1137 for c_vals in c: 

1138 tt, cc_val, kk = insert(x, [t, c_vals, k], m) 

1139 cc.append(cc_val) 

1140 return (tt, cc, kk) 

1141 else: 

1142 tt, cc, ier = _fitpack._insert(per, t, c, k, x, m) 

1143 if ier == 10: 

1144 raise ValueError("Invalid input data") 

1145 if ier: 

1146 raise TypeError("An error occurred") 

1147 return (tt, cc, k) 

1148 

1149 

1150def splder(tck, n=1): 

1151 """ 

1152 Compute the spline representation of the derivative of a given spline 

1153 

1154 Parameters 

1155 ---------- 

1156 tck : tuple of (t, c, k) 

1157 Spline whose derivative to compute 

1158 n : int, optional 

1159 Order of derivative to evaluate. Default: 1 

1160 

1161 Returns 

1162 ------- 

1163 tck_der : tuple of (t2, c2, k2) 

1164 Spline of order k2=k-n representing the derivative 

1165 of the input spline. 

1166 

1167 Notes 

1168 ----- 

1169 

1170 .. versionadded:: 0.13.0 

1171 

1172 See Also 

1173 -------- 

1174 splantider, splev, spalde 

1175 

1176 Examples 

1177 -------- 

1178 This can be used for finding maxima of a curve: 

1179 

1180 >>> from scipy.interpolate import splrep, splder, sproot 

1181 >>> x = np.linspace(0, 10, 70) 

1182 >>> y = np.sin(x) 

1183 >>> spl = splrep(x, y, k=4) 

1184 

1185 Now, differentiate the spline and find the zeros of the 

1186 derivative. (NB: `sproot` only works for order 3 splines, so we 

1187 fit an order 4 spline): 

1188 

1189 >>> dspl = splder(spl) 

1190 >>> sproot(dspl) / np.pi 

1191 array([ 0.50000001, 1.5 , 2.49999998]) 

1192 

1193 This agrees well with roots :math:`\\pi/2 + n\\pi` of 

1194 :math:`\\cos(x) = \\sin'(x)`. 

1195 

1196 """ 

1197 if n < 0: 

1198 return splantider(tck, -n) 

1199 

1200 t, c, k = tck 

1201 

1202 if n > k: 

1203 raise ValueError(("Order of derivative (n = %r) must be <= " 

1204 "order of spline (k = %r)") % (n, tck[2])) 

1205 

1206 # Extra axes for the trailing dims of the `c` array: 

1207 sh = (slice(None),) + ((None,)*len(c.shape[1:])) 

1208 

1209 with np.errstate(invalid='raise', divide='raise'): 

1210 try: 

1211 for j in range(n): 

1212 # See e.g. Schumaker, Spline Functions: Basic Theory, Chapter 5 

1213 

1214 # Compute the denominator in the differentiation formula. 

1215 # (and append traling dims, if necessary) 

1216 dt = t[k+1:-1] - t[1:-k-1] 

1217 dt = dt[sh] 

1218 # Compute the new coefficients 

1219 c = (c[1:-1-k] - c[:-2-k]) * k / dt 

1220 # Pad coefficient array to same size as knots (FITPACK 

1221 # convention) 

1222 c = np.r_[c, np.zeros((k,) + c.shape[1:])] 

1223 # Adjust knots 

1224 t = t[1:-1] 

1225 k -= 1 

1226 except FloatingPointError: 

1227 raise ValueError(("The spline has internal repeated knots " 

1228 "and is not differentiable %d times") % n) 

1229 

1230 return t, c, k 

1231 

1232 

1233def splantider(tck, n=1): 

1234 """ 

1235 Compute the spline for the antiderivative (integral) of a given spline. 

1236 

1237 Parameters 

1238 ---------- 

1239 tck : tuple of (t, c, k) 

1240 Spline whose antiderivative to compute 

1241 n : int, optional 

1242 Order of antiderivative to evaluate. Default: 1 

1243 

1244 Returns 

1245 ------- 

1246 tck_ader : tuple of (t2, c2, k2) 

1247 Spline of order k2=k+n representing the antiderivative of the input 

1248 spline. 

1249 

1250 See Also 

1251 -------- 

1252 splder, splev, spalde 

1253 

1254 Notes 

1255 ----- 

1256 The `splder` function is the inverse operation of this function. 

1257 Namely, ``splder(splantider(tck))`` is identical to `tck`, modulo 

1258 rounding error. 

1259 

1260 .. versionadded:: 0.13.0 

1261 

1262 Examples 

1263 -------- 

1264 >>> from scipy.interpolate import splrep, splder, splantider, splev 

1265 >>> x = np.linspace(0, np.pi/2, 70) 

1266 >>> y = 1 / np.sqrt(1 - 0.8*np.sin(x)**2) 

1267 >>> spl = splrep(x, y) 

1268 

1269 The derivative is the inverse operation of the antiderivative, 

1270 although some floating point error accumulates: 

1271 

1272 >>> splev(1.7, spl), splev(1.7, splder(splantider(spl))) 

1273 (array(2.1565429877197317), array(2.1565429877201865)) 

1274 

1275 Antiderivative can be used to evaluate definite integrals: 

1276 

1277 >>> ispl = splantider(spl) 

1278 >>> splev(np.pi/2, ispl) - splev(0, ispl) 

1279 2.2572053588768486 

1280 

1281 This is indeed an approximation to the complete elliptic integral 

1282 :math:`K(m) = \\int_0^{\\pi/2} [1 - m\\sin^2 x]^{-1/2} dx`: 

1283 

1284 >>> from scipy.special import ellipk 

1285 >>> ellipk(0.8) 

1286 2.2572053268208538 

1287 

1288 """ 

1289 if n < 0: 

1290 return splder(tck, -n) 

1291 

1292 t, c, k = tck 

1293 

1294 # Extra axes for the trailing dims of the `c` array: 

1295 sh = (slice(None),) + (None,)*len(c.shape[1:]) 

1296 

1297 for j in range(n): 

1298 # This is the inverse set of operations to splder. 

1299 

1300 # Compute the multiplier in the antiderivative formula. 

1301 dt = t[k+1:] - t[:-k-1] 

1302 dt = dt[sh] 

1303 # Compute the new coefficients 

1304 c = np.cumsum(c[:-k-1] * dt, axis=0) / (k + 1) 

1305 c = np.r_[np.zeros((1,) + c.shape[1:]), 

1306 c, 

1307 [c[-1]] * (k+2)] 

1308 # New knots 

1309 t = np.r_[t[0], t, t[-1]] 

1310 k += 1 

1311 

1312 return t, c, k