L iu :dZddlZddlZddlZddlmZmZmZm Z m Z ddl m Z dZ dZdZd Zd Zd Zdddd ddd dZdddd ddd dZdZdZGddZdZGddZdZedzedzddZdZdddd dddddZdddd dddddZddddd dddd d!Zy)"a Replicate FITPACK's logic for constructing smoothing spline functions and curves. Currently provides analogs of splrep and splprep python routines, i.e. curfit.f and parcur.f routines (the drivers are fpcurf.f and fppara.f, respectively) The Fortran sources are from https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/ .. [1] P. Dierckx, "Algorithms for smoothing data with periodic and parametric splines, Computer Graphics and Image Processing", 20 (1982) 171-184. :doi:`10.1016/0146-664X(82)90043-0`. .. [2] P. Dierckx, "Curve and surface fitting with splines", Monographs on Numerical Analysis, Oxford University Press, 1993. .. [3] P. Dierckx, "An algorithm for smoothing, differentiation and integration of experimental data using spline functions", Journal of Computational and Applied Mathematics, vol. I, no 3, p. 165 (1975). https://doi.org/10.1016/0771-050X(75)90034-0 N) _not_a_knotmake_interp_splineBSplinefpcheck _lsq_solve_qr)_dierckxgMbP?c|dz}t|||||\}}}tj|}t|||}t ||||} | j } tj | rttd| | fS)Nr) rnpascontiguousarrayr_compute_residualssumisnan ValueError_iermesg) xytkww2_cspl residualsfps f/mnt/ssd/data/python-lab/Trading/venv/lib/python3.12/site-packages/scipy/interpolate/_fitpack_repro.py_get_residualsr *s ABAq!Q*GAq! QA !Q C"2s1vq1I B xx|!%% b=c<||z dzjd}||zS)Nr raxis)r)rsplxrdeltas rrrCs'Qh]  Q  'E :r!ctj||||}tj||}tj|d||||df}|S)arAdd a new knot. (Approximately) replicate FITPACK's logic: 1. split the `x` array into knot intervals, ``t(j+k) <= x(i) <= t(j+k+1)`` 2. find the interval with the maximum sum of residuals 3. insert a new knot into the middle of that interval. NB: a new knot is in fact an `x` value at the middle of the interval. So *the knots are a subset of `x`*. This routine is an analog of https://github.com/scipy/scipy/blob/v1.11.4/scipy/interpolate/fitpack/fpcurf.f#L190-L215 (cf _split function) and https://github.com/scipy/scipy/blob/v1.11.4/scipy/interpolate/fitpack/fpknot.f N)r fpknotr searchsortedr_)rrrrnew_knotidx_tt_news radd_knotr.HsN"q!Q 2H OOAx (E EE!FU)Xqy0 1E Lr!ctj|t}tj|t}|tj|t}natj|t}|jdk7rt d|jd|dkj r t d|jdk(s|jdkDrt d |jd t|}|r+|jdk7rPt d |jd |d |jdk7rt d |jd |d |dddf}|jd|jdk7r&t d|jd|jd |jd|jdk7r&t d|jd|jd |jdk7s|dd|ddkj r t dtj|}|dkrt d|| t|}| t|}|||||||fS)zACommon input validations for generate_knots and make_splrep. dtypeNr w.ndim = z not implemented yet.rzWeights must be non-negativer z y.ndim = z not supported (must be 1 or 2.)z% != 2 not supported with parametric =.z% != 1 not supported with parametric =z"Weights is incompatible: w.shape =z != z Data is incompatible: x.shape = z and y.shape = z(Expect `x` to be an ordered 1D sequence.z"`s` must be non-negative. Got s = ) r asarrayfloat ones_likendimranyboolshapeoperatorindexminmax)rrrrsxbxe parametrics r_validate_inputsrD`s& 1E"A 1E"Ay LL% ( JJq & 66Q; {*?@A A E;;=;< <vv{affqjIAFF;&FGHHj!J 66Q; {*PJ>QRST T 66Q; {*PJ>QRST T agJwwqzQWWQZ>QWWKtAGG9ANOOwwqzQWWQZ!@AA z V z V aAq"b  r!rrArBrr@nestc #K|dk(r || tdt||}|yt|||||||tj|dk(\}}}}}}}t ||||||||Ed{y7w)a Generate knot vectors until the Least SQuares (LSQ) criterion is satified. Parameters ---------- x, y : array_like The data points defining the curve ``y = f(x)``. w : array_like, optional Weights. xb : float, optional The boundary of the approximation interval. If None (default), is set to ``x[0]``. xe : float, optional The boundary of the approximation interval. If None (default), is set to ``x[-1]``. k : int, optional The spline degree. Default is cubic, ``k = 3``. s : float, optional The smoothing factor. Default is ``s = 0``. nest : int, optional Stop when at least this many knots are placed. Yields ------ t : ndarray Knot vectors with an increasing number of knots. The generator is finite: it stops when the smoothing critetion is satisfied, or when then number of knots exceeds the maximum value: the user-provided `nest` or `x.size + k + 1` --- which is the knot vector for the interpolating spline. Examples -------- Generate some noisy data and fit a sequence of LSQ splines: >>> import numpy as np >>> import matplotlib.pyplot as plt >>> from scipy.interpolate import make_lsq_spline, generate_knots >>> rng = np.random.default_rng(12345) >>> x = np.linspace(-3, 3, 50) >>> y = np.exp(-x**2) + 0.1 * rng.standard_normal(size=50) >>> knots = list(generate_knots(x, y, s=1e-10)) >>> for t in knots[::3]: ... spl = make_lsq_spline(x, y, t) ... xs = xs = np.linspace(-3, 3, 201) ... plt.plot(xs, spl(xs), '-', label=f'n = {len(t)}', lw=3, alpha=0.7) >>> plt.plot(x, y, 'o', label='data') >>> plt.plot(xs, np.exp(-xs**2), '--') >>> plt.legend() Note that increasing the number of knots make the result follow the data more and more closely. Also note that a step of the generator may add multiple knots: >>> [len(t) for t in knots] [8, 9, 10, 12, 16, 24, 40, 48, 52, 54] Notes ----- The routine generates successive knots vectors of increasing length, starting from ``2*(k+1)`` to ``len(x) + k + 1``, trying to make knots more dense in the regions where the deviation of the LSQ spline from data is large. When the maximum number of knots, ``len(x) + k + 1`` is reached (this happens when ``s`` is small and ``nest`` is large), the generator stops, and the last output is the knots for the interpolation with the not-a-knot boundary condition. Knots are located at data sites, unless ``k`` is even and the number of knots is ``len(x) + k + 1``. In that case, the last output of the generator has internal knots at Greville sites, ``(x[1:] + x[:-1]) / 2``. .. versionadded:: 1.15.0 rNzs == 0 is interpolation onlyr rCrF)rrrDr r8_generate_knots_impl) rrrrArBrr@rGrs rgenerate_knotsrKsZ Av  q};< < 1 , 1aAr2"''!*/Aq!Q2r$AqA"qADQQQsA-A7/A50A7c #:K|tz}|j} |t| |zdzd|zdz}n#|d|dzzkrtd|dd|dzzdd|dzz} | |zdz} t j |g|dzz|g|dzzzt } | jd} d }d }t| D]}| |}t||| || \}}||z }t||ks|dkry| | k(rd}n=||z }||kDrt|z|z ndz}t|dzt||dzd}t|D]_}t|| ||} | jd} | | k\rt||} | y| |k\r| y||dz ksMt||| || \}}ayw) Nrr rE`nest` too small: nest = < 2*(k+1) = r3r0rgr)TOLsizer?rr r5r6r;ranger absintr>r.r)rrrrArBrr@rGaccmnminnmaxrnrfpoldrrfpmsnplusr&npl1js rrJrJs c'C A |1q519acAg& !QU) 9$-1Q3yPQRS S a!e9D q519D B41:ac *%8A  A B E1X.?&q!QQ7 2Av IO  9EBJE05 3ut|e+,qDaT5!8Q!78Eu ?AAq),A ADy1%Dy519}-aAqA> 1+ ?3.?b s FFFcdd}t|dzD]}||z|k7s ||||||zz z}|S)Ng?r )rR)rir^rresr@s rproddrbOsJ C 1Q3Z# q5A: AaD1QqS6M "C# Jr!c |jd}|||z dz ||z }|d|zz dz }tj|dz |dzft}t |dz D]G}||zdz}t |dzD],}||z} || |zdz|| z t || ||z |||f<.I|||z |zz}tj t |dz D cgc]} | c} tj} ||z dz } || | fScc} w)aDiscontinuity matrix: jumps of k-th derivatives of b-splines at internal knots. See Eqs. (9)-(10) of Ref. [1], or, equivalently, Eq. (3.43) of Ref. [2]. This routine assumes internal knots are all simple (have multiplicity =1). Parameters ---------- t : ndarray, 1D, shape(n,) Knots. k : int The spline degree Returns ------- disc : ndarray, shape(n-2*k-1, k+2) The jumps of the k-th derivatives of b-splines at internal knots, ``t[k+1], ...., t[n-k-1]``. offset : ndarray, shape(2-2*k-1,) Offsets nc : int Notes ----- The normalization here follows FITPACK: (https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fpdisc.f#L36) The k-th derivative jumps are multiplied by a factor:: (delta / nrint)**k / k! where ``delta`` is the length of the interval spanned by internal knots, and ``nrint`` is one less the number of internal knots (i.e., the number of subintervals between them). References ---------- .. [1] Paul Dierckx, Algorithms for smoothing data with periodic and parametric splines, Computer Graphics and Image Processing, vol. 20, p. 171 (1982). :doi:`10.1016/0146-664X(82)90043-0` .. [2] Tom Lyche and Knut Morken, Spline methods, http://www.uio.no/studier/emner/matnat/ifi/INF-MAT5340/v05/undervisningsmateriale/ rrr r0)r;r emptyr6rRrbarrayint64) rrrYr&nrintmatrjjr^iir`offsetncs rdiscrmWs4^  A a!eaiL1Q4 E !GaKE 88UQYA&e 4DEAIE FQJA, EBRAa!eaiL1Q4/5Aq!3DDDRL EE UE\A DXX%a.1Qq1 BF QB  2s C>c&eZdZdZdddddZdZy)Fa The r.h.s. of ``f(p) = s``. Given scalar `p`, we solve the system of equations in the LSQ sense: | A | @ | c | = | y | | B / p | | 0 | | 0 | where `A` is the matrix of b-splines and `b` is the discontinuity matrix (the jumps of the k-th derivatives of b-spline basis elements at knots). Since we do that repeatedly while minimizing over `p`, we QR-factorize `A` only once and update the QR factorization only of the `B` rows of the augmented matrix |A, B/p|. The system of equations is Eq. (15) Ref. [1]_, the strategy and implementation follows that of FITPACK, see specific links below. References ---------- [1] P. Dierckx, Algorithms for Smoothing Data with Periodic and Parametric Splines, COMPUTER GRAPHICS AND IMAGE PROCESSING vol. 20, pp 171-184 (1982.) https://doi.org/10.1016/0146-664X(82)90043-0 N)RYc||_||_||_||_|t j |t n|}|jdk7rtd|jd||_ ||_ |jdk7rtd|jdt||\} } } ||t|||||\}}} |jd|z dz } |dz}|jd|k7r"td |jdd |dzd t j| jd|jdft }tj|d| |f|_t j| | jdz|jdzft }|d| ddf|d| d|f<||_tjt j$| tj&| f|_| |_| |_y) Nr0rr2z != 1.r z&F: expected y.ndim == 2, got y.ndim = z instead.rzInternal error: R.shape[1] =z != k+1 =r3)rrrrr r7r6r8rrr@rmrr;zerosr*YYAAarangerfrkrlb)selfrrrrr@rrprqrwb_offsetb_ncrrlnzzrus r__init__z F.__init__s,-IBLL% (1 66Q; {&12 2 66Q;FQVVKyQR R!AJ8T 9#Aq!Q2GAq! WWQZ!^a  U 771: < ~Z!A#JK K HHaggaj!''!*-U ;%%#2 "XXrAGGAJq1 ?"ay3B38 eeBIIb98CD r!c>|jj}|jj}|j}|j|z ||dddf<|j j}t j|||||t j|||}t|j||j}t|jdz||j|j}|j!} ||_| |j$z S)N)startrowr )rucopyrkrlrwrtr qr_reducefpbackrrrrrrrrrr@) rxpABrkrlQYrrrrs r__call__z F.__call__s WW\\^!!# WWVVaZ236 WW\\^ 2vr2; OOBB 'dffa(&tvvqy#dff+tvvF ]]_DFF{r!)N)__name__ __module__ __qualname____doc__r}rr!rroros0/44/br!roc|||z z}|||z z}|||z z}|tjk(r||z||zz |z S||z|z||z|zz||z|zz ||z||zz||zzz S)aThe root of r(p) = (u*p + v) / (p + w) given three points and values, (p1, f2), (p2, f2) and (p3, f3). The FITPACK analog adjusts the bounds, and we do not https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fprati.f NB: FITPACK uses p < 0 to encode p=infinity. We just use the infinity itself. Since the bracket is ``p1 <= p2 <= p3``, ``p3`` can be infinite (in fact, this is what the minimizer starts with, ``p3=inf``). )r inf) p1f1p2f2p3f3h1h2h3s rfpratirs rBwB rBwB rBwB RVV|BB"$$ U2X2b 2b58 + ,22 20E FFr!ceZdZdZy)Bunchc <|jjdi|y)Nr)__dict__update)rxkwargss rr}zBunch.__init__s &v&r!N)rrrr}rr!rrrs'r!rzerror. a theoretically impossible result was found during the iteration process for finding a smoothing spline with fp = s. probably causes : s too small. z6the weighted sum of squared residuals is becoming NaN zthere is an approximation returned but the corresponding weighted sum of squared residuals does not satisfy the condition abs(fp-s)/s < tol. a7error. the maximal number of iterations maxit (set to 20 by the program) allowed for finding a smoothing spline with fp=s has been reached. probably causes : s too small there is an approximation returned but the corresponding weighted sum of squared residuals does not satisfy the condition abs(fp-s)/s < tol. )rr rEc Td}d}d}d\}}|\\} } \} } |} ttD]}| || }}t||krd\}}n|dk(r)|| z |kr|} |} | |z} | | kr | |z||zz} I|dkrd}|dk(r<| |z |kr-|} |} | |z } | tjk7r| | kr ||z| |zz} |dkDrd}| |ks|| krd\}}n%t | | ||| | } |dkr||} } ||} } d \}}|dk7r't jtt|d t|| | S) aSolve `f(p) = 0` using a rational function approximation. In a nutshell, since the function f(p) is known to be monotonically decreasing, we - maintain the bracket (p1, f1), (p2, f2) and (p3, f3) - at each iteration step, approximate f(p) by a rational function r(p) = (u*p + v) / (p + w) and make a step to p_new to the root of f(p): r(p_new) = 0. The coefficients u, v and w are found from the bracket values p1..3 and f1...3 The algorithm and implementation follows https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fpcurf.f#L229 and https://github.com/scipy/scipy/blob/maintenance/1.11.x/scipy/interpolate/fitpack/fppara.f#L290 Note that the latter is for parametric splines and the former is for 1D spline functions. The minimization is indentical though [modulo a summation over the dimensions in the computation of f(p)], so we reuse the minimizer for both d=1 and d>1. g?g?g{Gz?)rr)rTrr)r F)rEFr ) stacklevel) convergedroot iterationsier) rRMAXITrSr rrwarningswarnRuntimeWarningrr)fp0bracketrUcon1con9con4ich1ich3rrrrritrrrrs r root_ratir4s, D D DJD$!HRhr2 AEl4"AaDB r7S=$NC  19Bw#~dF7D2d7*Q6D 19Bw#~dF 0` is the input parameter. In other words, we balance maximizing the smoothness (measured as the jumps of the derivative, the first criterion), and the deviation of :math:`g(x_j)` from the data :math:`y_j` (the second criterion). Note that the summation in the second criterion is over all data points, and in the first criterion it is over the internal spline knots (i.e. those with ``xb < t[i] < xe``). The spline knots are in general a subset of data, see `generate_knots` for details. Also note the difference of this routine to `make_lsq_spline`: the latter routine does not consider smoothness and simply solves a least-squares problem .. math:: \sum w_j \times (g(x_j) - y_j)^2 \to \mathrm{min} for a spline function :math:`g(x)` with a _fixed_ knot vector ``t``. .. versionadded:: 1.15.0 rNs==0 is for interpolation only)rFrIr)rrrDrr) rrrrArBrr@rrGrs r make_splreprsn Av =AMT-==> >!!Q!,,,Q1aBuUAq!Q2r AqA"qA NC EE!Q$KCE Jr!)ruubuerr@rrGc >tj|d}|h|ddddf|ddddfz dz} tj| jdj }tj d||dz f}|dk(r,||| t dt||j|d|fSt|||||||d \}}}}}}}t||||||||| } | jj} t| j| | jd} | |fS) ab Create a smoothing parametric B-spline curve with bounded error, minimizing derivative jumps. Given a list of N 1D arrays, `x`, which represent a curve in N-dimensional space parametrized by `u`, find a smooth approximating spline curve ``g(u)``. Parameters ---------- x : array_like, shape (m, ndim) Sampled data points representing the curve in ``ndim`` dimensions. The typical use is a list of 1D arrays, each of length ``m``. w : array_like, shape(m,), optional Strictly positive 1D array of weights. The weights are used in computing the weighted least-squares spline fit. If the errors in the `x` values have standard deviation given by the vector d, then `w` should be 1/d. Default is ``np.ones(m)``. u : array_like, optional An array of parameter values for the curve in the parametric form. If not given, these values are calculated automatically, according to:: v[0] = 0 v[i] = v[i-1] + distance(x[i], x[i-1]) u[i] = v[i] / v[-1] ub, ue : float, optional The end-points of the parameters interval. Default to ``u[0]`` and ``u[-1]``. k : int, optional Degree of the spline. Cubic splines, ``k=3``, are recommended. Even values of `k` should be avoided especially with a small ``s`` value. Default is ``k=3`` s : float, optional A smoothing condition. The amount of smoothness is determined by satisfying the conditions:: sum((w * (g(u) - x))**2) <= s, where ``g(u)`` is the smoothed approximation to ``x``. The user can use `s` to control the trade-off between closeness and smoothness of fit. Larger ``s`` means more smoothing while smaller values of ``s`` indicate less smoothing. Recommended values of ``s`` depend on the weights, ``w``. If the weights represent the inverse of the standard deviation of ``x``, then a good ``s`` value should be found in the range ``(m - sqrt(2*m), m + sqrt(2*m))``, where ``m`` is the number of data points in ``x`` and ``w``. t : array_like, optional The spline knots. If None (default), the knots will be constructed automatically. There must be at least ``2*k + 2`` and at most ``m + k + 1`` knots. nest : int, optional The target length of the knot vector. Should be between ``2*(k + 1)`` (the minimum number of knots for a degree-``k`` spline), and ``m + k + 1`` (the number of knots of the interpolating spline). The actual number of knots returned by this routine may be slightly larger than `nest`. Default is None (no limit, add up to ``m + k + 1`` knots). Returns ------- spl : a `BSpline` instance For `s=0`, ``spl(u) == x``. For non-zero values of ``s``, `spl` represents the smoothed approximation to ``x``, generally with fewer knots. u : ndarray The values of the parameters See Also -------- generate_knots : is used under the hood for generating the knots make_splrep : the analog of this routine 1D functions make_interp_spline : construct an interpolating spline (``s = 0``) make_lsq_spline : construct the least-squares spline given the knot vector splprep : a FITPACK analog of this routine Notes ----- Given a set of :math:`m` data points in :math:`D` dimensions, :math:`\vec{x}_j`, with :math:`j=1, ..., m` and :math:`\vec{x}_j = (x_{j; 1}, ..., x_{j; D})`, this routine constructs the parametric spline curve :math:`g_a(u)` with :math:`a=1, ..., D`, to minimize the sum of jumps, :math:`D_{i; a}`, of the ``k``-th derivative at the internal knots (:math:`u_b < t_i < u_e`), where .. math:: D_{i; a} = g_a^{(k)}(t_i + 0) - g_a^{(k)}(t_i - 0) Specifically, the routine constructs the spline function :math:`g(u)` which minimizes .. math:: \sum_i \sum_{a=1}^D | D_{i; a} |^2 \to \mathrm{min} provided that .. math:: \sum_{j=1}^m \sum_{a=1}^D (w_j \times (g_a(u_j) - x_{j; a}))^2 \leqslant s where :math:`u_j` is the value of the parameter corresponding to the data point :math:`(x_{j; 1}, ..., x_{j; D})`, and :math:`s > 0` is the input parameter. In other words, we balance maximizing the smoothness (measured as the jumps of the derivative, the first criterion), and the deviation of :math:`g(u_j)` from the data :math:`x_j` (the second criterion). Note that the summation in the second criterion is over all data points, and in the first criterion it is over the internal spline knots (i.e. those with ``ub < t[i] < ue``). The spline knots are in general a subset of data, see `generate_knots` for details. .. versionadded:: 1.15.0 References ---------- .. [1] P. Dierckx, "Algorithms for smoothing data with periodic and parametric splines, Computer Graphics and Image Processing", 20 (1982) 171-184. .. [2] P. Dierckx, "Curve and surface fitting with splines", Monographs on Numerical Analysis, Oxford University Press, 1993. rr#Nr4r rr)rr$TrIr)r stacksqrtrcumsumr*rrTrDrrrrr) rrrrrrr@rrGdprccspl1s r make_splpreprTs)t A yAh3B36"Q & GGRHH!H$ % , , . EE!Q2Y, Av =AMT-==> >!!QSSAA699,Q1aBtTAq!Q2r AqA"qA NC B 355"cee! ,D 7Nr!) rrr<numpyr _bsplinesrrrrrr rPrr rr.rDrKrJrbrmrorr _iermesg1rrrrrrr!rrs&  2 0.!b#tQTXRv%)Tda14I FFRccLG&''      "ZFz"&$41TPT;| DTQ!t$BJ$41TPTPr!