-
Notifications
You must be signed in to change notification settings - Fork 26
/
_sspumv.py
369 lines (275 loc) · 12.2 KB
/
_sspumv.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
# -*- coding: utf-8 -*-
"""
Univariate/multivariate cubic smoothing spline implementation
"""
import functools
from typing import Optional, Union, Tuple, List
import numpy as np
import scipy.sparse as sp
import scipy.sparse.linalg as la
from scipy.interpolate import PPoly
from ._base import ISplinePPForm, ISmoothingSpline
from ._types import UnivariateDataType, MultivariateDataType
from ._reshape import to_2d, prod
class SplinePPForm(ISplinePPForm[np.ndarray, int], PPoly):
"""The base class for univariate/multivariate spline in piecewise polynomial form
Piecewise polynomial in terms of coefficients and breakpoints.
Notes
-----
Inherited from :py:class:`scipy.interpolate.PPoly`
"""
__module__ = 'csaps'
@property
def breaks(self) -> np.ndarray:
return self.x
@property
def coeffs(self) -> np.ndarray:
return self.c
@property
def gcv(self) -> float:
return self.gcv
@property
def order(self) -> int:
return self.c.shape[0]
@property
def pieces(self) -> int:
return self.c.shape[1]
@property
def ndim(self) -> int:
"""Returns the number of spline dimensions
The number of dimensions is product of shape without ``shape[self.axis]``.
"""
shape = list(self.shape)
shape.pop(self.axis)
return prod(shape)
@property
def shape(self) -> Tuple[int]:
"""Returns the source data shape
"""
shape: List[int] = list(self.c.shape[2:])
shape.insert(self.axis, self.c.shape[1] + 1)
return tuple(shape)
def __repr__(self): # pragma: no cover
return (
f'{type(self).__name__}\n'
f' breaks: {self.breaks}\n'
f' coeffs shape: {self.coeffs.shape}\n'
f' data shape: {self.shape}\n'
f' axis: {self.axis}\n'
f' pieces: {self.pieces}\n'
f' order: {self.order}\n'
f' ndim: {self.ndim}\n'
f' gcv: {self.gcv}\n'
)
class CubicSmoothingSpline(ISmoothingSpline[
SplinePPForm,
float,
UnivariateDataType,
int,
Union[bool, str]
]):
"""Cubic smoothing spline
The cubic spline implementation for univariate/multivariate data.
Parameters
----------
xdata : np.ndarray, sequence, vector-like
X input 1-D data vector (data sites: ``x1 < x2 < ... < xN``)
ydata : np.ndarray, vector-like, sequence[vector-like]
Y input 1-D data vector or ND-array with shape[axis] equal of `xdata` size)
weights : [*Optional*] np.ndarray, list
Weights 1-D vector with size equal of ``xdata`` size
smooth : [*Optional*] float
Smoothing parameter in range [0, 1] where:
- 0: The smoothing spline is the least-squares straight line fit
- 1: The cubic spline interpolant with natural condition
axis : [*Optional*] int
Axis along which ``ydata`` is assumed to be varying.
Meaning that for x[i] the corresponding values are np.take(ydata, i, axis=axis).
By default is -1 (the last axis).
normalizedsmooth : [*Optional*] bool
If True, the smooth parameter is normalized such that results are invariant to xdata range
and less sensitive to nonuniformity of weights and xdata clumping
.. versionadded:: 1.1.0
calculate_gcv : [*Optional*] bool
If True, the Generalized Cross Validation criterion value will be calculated for the spline upon creation.
The GCV can be minimized to attempt to identify an optimal smoothing parameter, while penalizing over fitting.
Strictly speaking this involves generating splines for all smoothing parameter values across the valid domain
of the smoothing parameter [0,1] then selecting the smoothing parameter that produces the lowest GCV.
[See GCV in TEOSL (page 244 section 7.52)](https://hastie.su.domains/ElemStatLearn/printings/ESLII_print12_toc.pdf) for more information on methodology.
The simplest way to use the GCV is to setup a loop that generates a spline with your data at every step, a step size
of 0.001 is generally sufficient, and keep track of the spline that produces the lowest GCV. Setting this parameter to True,
does _not_ generate the lowest GCV, it simply generates the value for the spline with this specific smoothing parameter,
so that you may decide how to optimize it for yourself.
"""
__module__ = 'csaps'
def __init__(self,
xdata: UnivariateDataType,
ydata: MultivariateDataType,
weights: Optional[UnivariateDataType] = None,
smooth: Optional[float] = None,
axis: int = -1,
normalizedsmooth: bool = False,
calculate_gcv: bool = False):
x, y, w, shape, axis = self._prepare_data(xdata, ydata, weights, axis)
if calculate_gcv:
coeffs, smooth, gcv = self._make_spline(x, y, w, smooth, shape, normalizedsmooth, calculate_gcv)
self.gcv = gcv
else:
coeffs, smooth = self._make_spline(x, y, w, smooth, shape, normalizedsmooth, calculate_gcv)
self.gcv = None
spline = SplinePPForm.construct_fast(coeffs, x, axis=axis)
self._smooth = smooth
self._spline = spline
def __call__(self,
x: UnivariateDataType,
nu: Optional[int] = None,
extrapolate: Optional[Union[bool, str]] = None) -> np.ndarray:
"""Evaluate the spline for given data
Parameters
----------
x : 1-d array-like
Points to evaluate the spline at.
nu : [*Optional*] int
Order of derivative to evaluate. Must be non-negative.
extrapolate : [*Optional*] bool or 'periodic'
If bool, determines whether to extrapolate to out-of-bounds points
based on first and last intervals, or to return NaNs. If 'periodic',
periodic extrapolation is used. Default is True.
Notes
-----
Derivatives are evaluated piecewise for each polynomial
segment, even if the polynomial is not differentiable at the
breakpoints. The polynomial intervals are considered half-open,
``[a, b)``, except for the last interval which is closed
``[a, b]``.
"""
if nu is None:
nu = 0
return self._spline(x, nu=nu, extrapolate=extrapolate)
@property
def smooth(self) -> float:
"""Returns the smoothing factor
Returns
-------
smooth : float
Smoothing factor in the range [0, 1]
"""
return self._smooth
@property
def spline(self) -> SplinePPForm:
"""Returns the spline description in `SplinePPForm` instance
Returns
-------
spline : SplinePPForm
The spline representation in :class:`SplinePPForm` instance
"""
return self._spline
@staticmethod
def _prepare_data(xdata, ydata, weights, axis):
xdata = np.asarray(xdata, dtype=np.float64)
ydata = np.asarray(ydata, dtype=np.float64)
if xdata.ndim > 1:
raise ValueError("'xdata' must be a vector")
if xdata.size < 2:
raise ValueError("'xdata' must contain at least 2 data points.")
axis = ydata.ndim + axis if axis < 0 else axis
if ydata.shape[axis] != xdata.size:
raise ValueError(
f"'ydata' data must be a 1-D or N-D array with shape[{axis}] "
f"that is equal to 'xdata' size ({xdata.size})")
# Rolling axis for using its shape while constructing coeffs array
shape = np.rollaxis(ydata, axis).shape
# Reshape ydata N-D array to 2-D NxM array where N is the data
# dimension and M is the number of data points.
ydata = to_2d(ydata, axis)
if weights is None:
weights = np.ones_like(xdata)
else:
weights = np.asarray(weights, dtype=np.float64)
if weights.size != xdata.size:
raise ValueError('Weights vector size must be equal of xdata size')
return xdata, ydata, weights, shape, axis
@staticmethod
def _compute_smooth(a, b):
"""
The calculation of the smoothing spline requires the solution of a
linear system whose coefficient matrix has the form p*A + (1-p)*B, with
the matrices A and B depending on the data sites x. The default value
of p makes p*trace(A) equal (1 - p)*trace(B).
"""
def trace(m: sp.dia_matrix):
return m.diagonal().sum()
return 1. / (1. + trace(a) / (6. * trace(b)))
@staticmethod
def _normalize_smooth(x: np.ndarray, w: np.ndarray, smooth: Optional[float]):
"""
See the explanation here: https://github.com/espdev/csaps/pull/47
"""
span = np.ptp(x)
eff_x = 1 + (span ** 2) / np.sum(np.diff(x) ** 2)
eff_w = np.sum(w) ** 2 / np.sum(w ** 2)
k = 80 * (span ** 3) * (x.size ** -2) * (eff_x ** -0.5) * (eff_w ** -0.5)
s = 0.5 if smooth is None else smooth
p = s / (s + (1 - s) * k)
return p
@staticmethod
def _make_spline(x, y, w, smooth, shape, normalizedsmooth, calculate_gcv):
pcount = x.size
dx = np.diff(x)
if not all(dx > 0): # pragma: no cover
raise ValueError(
"Items of 'xdata' vector must satisfy the condition: x1 < x2 < ... < xN")
dy = np.diff(y, axis=1)
dy_dx = dy / dx
if pcount == 2:
# The corner case for the data with 2 points (1 breaks interval)
# In this case we have 2-ordered spline and linear interpolation in fact
yi = y[:, 0][:, np.newaxis]
c_shape = (2, pcount - 1) + shape[1:]
c = np.vstack((dy_dx, yi)).reshape(c_shape)
p = 1.0
return c, p
# Create diagonal sparse matrices
diags_r = np.vstack((dx[1:], 2 * (dx[1:] + dx[:-1]), dx[:-1]))
r = sp.spdiags(diags_r, [-1, 0, 1], pcount - 2, pcount - 2)
dx_recip = 1. / dx
diags_qtw = np.vstack((dx_recip[:-1], -(dx_recip[1:] + dx_recip[:-1]), dx_recip[1:]))
diags_sqrw_recip = 1. / np.sqrt(w)
# Calculate and store qt separately, so we can use it later to calculate the degrees of freedom for the gcv
qt = sp.diags(diags_qtw, [0, 1, 2], (pcount - 2, pcount))
qtw = ( qt @ sp.diags(diags_sqrw_recip, 0, (pcount, pcount)))
qtw = qtw @ qtw.T
p = smooth
if normalizedsmooth:
p = CubicSmoothingSpline._normalize_smooth(x, w, smooth)
elif smooth is None:
p = CubicSmoothingSpline._compute_smooth(r, qtw)
pp = (6. * (1. - p))
# Solve linear system for the 2nd derivatives
a = pp * qtw + p * r
b = np.diff(dy_dx, axis=1).T
u = la.spsolve(a, b)
if u.ndim < 2:
u = u[np.newaxis]
if y.shape[0] == 1:
u = u.T
dx = dx[:, np.newaxis]
vpad = functools.partial(np.pad, pad_width=[(1, 1), (0, 0)], mode='constant')
d1 = np.diff(vpad(u), axis=0) / dx
d2 = np.diff(vpad(d1), axis=0)
diags_w_recip = 1. / w
w = sp.diags(diags_w_recip, 0, (pcount, pcount))
yi = y.T - (pp * w) @ d2
pu = vpad(p * u)
c1 = np.diff(pu, axis=0) / dx
c2 = 3. * pu[:-1, :]
c3 = np.diff(yi, axis=0) / dx - dx * (2. * pu[:-1, :] + pu[1:, :])
c4 = yi[:-1, :]
c_shape = (4, pcount - 1) + shape[1:]
c = np.vstack((c1, c2, c3, c4)).reshape(c_shape)
# As calculating the GCV is a relatively expensive operation, it's best to add the calculation inline and make it optional
if calculate_gcv:
degrees_of_freedom = p * (la.inv(p * W + ((1-p) * 6 * qt.T @ la.inv(r) @ qt)) @ W).trace()
gcv = np.linalg.norm(y-y_pred)**2 / (y.size - degrees_of_freedom)**2
return c, p, gcv
return c, p