/
_splines.py
384 lines (304 loc) · 11.9 KB
/
_splines.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
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
import cupy
from cupy._core._scalar import get_typename
from cupy_backends.cuda.api import runtime
from cupyx.scipy.signal._signaltools import lfilter
if runtime.is_hip:
SYMIIR2_KERNEL = r"""#include <hip/hip_runtime.h>
"""
else:
SYMIIR2_KERNEL = r"""
#include <cuda_runtime.h>
#include <device_launch_parameters.h>
"""
SYMIIR2_KERNEL = SYMIIR2_KERNEL + r"""
#include <cupy/math_constants.h>
#include <cupy/carray.cuh>
template<typename T>
__device__ T _compute_symiirorder2_fwd_hc(
const int k, const T cs, const T r, const T omega) {
T base;
if(k < 0) {
return 0;
}
if(omega == 0.0) {
base = cs * pow(r, ((T) k)) * (k + 1);
} else if(omega == M_PI) {
base = cs * pow(r, ((T) k)) * (k + 1) * (1 - 2 * (k % 2));
} else {
base = (cs * pow(r, ((T) k)) * sin(omega * (k + 1)) /
sin(omega));
}
return base;
}
template<typename T>
__global__ void compute_symiirorder2_fwd_sc(
const int n, const int off, const T* cs_ptr, const T* r_ptr,
const T* omega_ptr, const double precision, bool* valid, T* out) {
int idx = blockDim.x * blockIdx.x + threadIdx.x;
if(idx + off >= n) {
return;
}
const T cs = cs_ptr[0];
const T r = r_ptr[0];
const T omega = omega_ptr[0];
T val = _compute_symiirorder2_fwd_hc<T>(idx + off + 1, cs, r, omega);
T err = val * val;
out[idx] = val;
valid[idx] = err <= precision;
}
template<typename T>
__device__ T _compute_symiirorder2_bwd_hs(
const int ki, const T cs, const T rsq, const T omega) {
T c0;
T gamma;
T cssq = cs * cs;
int k = abs(ki);
T rsupk = pow(rsq, ((T) k) / ((T) 2.0));
if(omega == 0.0) {
c0 = (1 + rsq) / ((1 - rsq) * (1 - rsq) * (1 - rsq)) * cssq;
gamma = (1 - rsq) / (1 + rsq);
return c0 * rsupk * (1 + gamma * k);
}
if(omega == M_PI) {
c0 = (1 + rsq) / ((1 - rsq) * (1 - rsq) * (1 - rsq)) * cssq;
gamma = (1 - rsq) / (1 + rsq) * (1 - 2 * (k % 2));
return c0 * rsupk * (1 + gamma * k);
}
c0 = (cssq * (1.0 + rsq) / (1.0 - rsq) /
(1 - 2 * rsq * cos(2 * omega) + rsq * rsq));
gamma = (1.0 - rsq) / (1.0 + rsq) / tan(omega);
return c0 * rsupk * (cos(omega * k) + gamma * sin(omega * k));
}
template<typename T>
__global__ void compute_symiirorder2_bwd_sc(
const int n, const int off, const int l_off, const int r_off,
const T* cs_ptr, const T* rsq_ptr, const T* omega_ptr,
const double precision, bool* valid, T* out) {
int idx = blockDim.x * blockIdx.x + threadIdx.x;
if(idx + off >= n) {
return;
}
const T cs = cs_ptr[0];
const T rsq = rsq_ptr[0];
const T omega = omega_ptr[0];
T v1 = _compute_symiirorder2_bwd_hs<T>(idx + l_off + off, cs, rsq, omega);
T v2 = _compute_symiirorder2_bwd_hs<T>(idx + r_off + off, cs, rsq, omega);
T diff = v1 + v2;
T err = diff * diff;
out[idx] = diff;
valid[idx] = err <= precision;
}
"""
SYMIIR2_MODULE = cupy.RawModule(
code=SYMIIR2_KERNEL, options=('-std=c++11',),
name_expressions=[f'compute_symiirorder2_bwd_sc<{t}>'
for t in ['float', 'double']] +
[f'compute_symiirorder2_fwd_sc<{t}>'
for t in ['float', 'double']])
def _get_module_func(module, func_name, *template_args):
args_dtypes = [get_typename(arg.dtype) for arg in template_args]
template = ', '.join(args_dtypes)
kernel_name = f'{func_name}<{template}>' if template_args else func_name
kernel = module.get_function(kernel_name)
return kernel
def _find_initial_cond(all_valid, cum_poly, n, off=0):
indices = cupy.where(all_valid)[0] + 1 + off
zi = cupy.nan
if indices.size > 0:
zi = cupy.where(
indices[0] >= n, cupy.nan, cum_poly[indices[0] - 1 - off])
return zi
def symiirorder1(input, c0, z1, precision=-1.0):
"""
Implement a smoothing IIR filter with mirror-symmetric boundary conditions
using a cascade of first-order sections. The second section uses a
reversed sequence. This implements a system with the following
transfer function and mirror-symmetric boundary conditions::
c0
H(z) = ---------------------
(1-z1/z) (1 - z1 z)
The resulting signal will have mirror symmetric boundary conditions
as well.
Parameters
----------
input : ndarray
The input signal.
c0, z1 : scalar
Parameters in the transfer function.
precision :
Specifies the precision for calculating initial conditions
of the recursive filter based on mirror-symmetric input.
Returns
-------
output : ndarray
The filtered signal.
"""
if cupy.abs(z1) >= 1:
raise ValueError('|z1| must be less than 1.0')
if precision <= 0.0 or precision > 1.0:
precision = cupy.finfo(input.dtype).resolution
precision *= precision
pos = cupy.arange(1, input.size + 1, dtype=input.dtype)
pow_z1 = z1 ** pos
diff = pow_z1 * cupy.conjugate(pow_z1)
cum_poly = cupy.cumsum(pow_z1 * input) + input[0]
all_valid = diff <= precision
zi = _find_initial_cond(all_valid, cum_poly, input.size)
if cupy.isnan(zi):
raise ValueError(
'Sum to find symmetric boundary conditions did not converge.')
# Apply first the system 1 / (1 - z1 * z^-1)
y1, _ = lfilter(
cupy.ones(1, dtype=input.dtype), cupy.r_[1, -z1], input[1:], zi=zi)
y1 = cupy.r_[zi, y1]
# Compute backward symmetric condition and apply the system
# c0 / (1 - z1 * z)
zi = -c0 / (z1 - 1.0) * y1[-1]
out, _ = lfilter(c0, cupy.r_[1, -z1], y1[:-1][::-1], zi=zi)
return cupy.r_[out[::-1], zi]
def _compute_symiirorder2_fwd_hc(k, cs, r, omega):
base = None
if omega == 0.0:
base = cs * cupy.power(r, k) * (k+1)
elif omega == cupy.pi:
base = cs * cupy.power(r, k) * (k + 1) * (1 - 2 * (k % 2))
else:
base = (cs * cupy.power(r, k) * cupy.sin(omega * (k + 1)) /
cupy.sin(omega))
return cupy.where(k < 0, 0.0, base)
def _compute_symiirorder2_bwd_hs(k, cs, rsq, omega):
cssq = cs * cs
k = cupy.abs(k)
rsupk = cupy.power(rsq, k / 2.0)
if omega == 0.0:
c0 = (1 + rsq) / ((1 - rsq) * (1 - rsq) * (1 - rsq)) * cssq
gamma = (1 - rsq) / (1 + rsq)
return c0 * rsupk * (1 + gamma * k)
if omega == cupy.pi:
c0 = (1 + rsq) / ((1 - rsq) * (1 - rsq) * (1 - rsq)) * cssq
gamma = (1 - rsq) / (1 + rsq) * (1 - 2 * (k % 2))
return c0 * rsupk * (1 + gamma * k)
c0 = (cssq * (1.0 + rsq) / (1.0 - rsq) /
(1 - 2 * rsq * cupy.cos(2 * omega) + rsq * rsq))
gamma = (1.0 - rsq) / (1.0 + rsq) / cupy.tan(omega)
return c0 * rsupk * (cupy.cos(omega * k) + gamma * cupy.sin(omega * k))
def symiirorder2(input, r, omega, precision=-1.0):
"""
Implement a smoothing IIR filter with mirror-symmetric boundary conditions
using a cascade of second-order sections. The second section uses a
reversed sequence. This implements the following transfer function::
cs^2
H(z) = ---------------------------------------
(1 - a2/z - a3/z^2) (1 - a2 z - a3 z^2 )
where::
a2 = 2 * r * cos(omega)
a3 = - r ** 2
cs = 1 - 2 * r * cos(omega) + r ** 2
Parameters
----------
input : ndarray
The input signal.
r, omega : float
Parameters in the transfer function.
precision : float
Specifies the precision for calculating initial conditions
of the recursive filter based on mirror-symmetric input.
Returns
-------
output : ndarray
The filtered signal.
"""
if r >= 1.0:
raise ValueError('r must be less than 1.0')
if precision <= 0.0 or precision > 1.0:
if input.dtype is cupy.dtype(cupy.float64):
precision = 1e-11
elif input.dtype is cupy.dtype(cupy.float32):
precision = 1e-6
else:
precision = 10 ** -cupy.finfo(input.dtype).iexp
block_sz = 128
rsq = r * r
a2 = 2 * r * cupy.cos(omega)
a3 = -rsq
cs = cupy.atleast_1d(1 - 2 * r * cupy.cos(omega) + rsq)
omega = cupy.asarray(omega, cs.dtype)
r = cupy.asarray(r, cs.dtype)
rsq = cupy.asarray(rsq, cs.dtype)
precision *= precision
# First compute the symmetric forward starting conditions
compute_symiirorder2_fwd_sc = _get_module_func(
SYMIIR2_MODULE, 'compute_symiirorder2_fwd_sc', cs)
diff = cupy.empty((block_sz + 1,), dtype=cs.dtype)
all_valid = cupy.empty((block_sz + 1,), dtype=cupy.bool_)
starting_diff = cupy.arange(2, dtype=input.dtype)
starting_diff = _compute_symiirorder2_fwd_hc(starting_diff, cs, r, omega)
y0 = cupy.nan
y1 = cupy.nan
for i in range(0, input.size + 2, block_sz):
compute_symiirorder2_fwd_sc(
(1,), (block_sz + 1,), (
input.size + 2, i, cs, r, omega, precision, all_valid, diff))
input_slice = input[i:i + block_sz]
diff_y0 = diff[:-1][:input_slice.size]
diff_y1 = diff[1:][:input_slice.size]
if cupy.isnan(y0):
cum_poly_y0 = cupy.cumsum(
diff_y0 * input_slice) + starting_diff[0] * input[0]
y0 = _find_initial_cond(
all_valid[:-1][:input_slice.size], cum_poly_y0, input.size, i)
if cupy.isnan(y1):
cum_poly_y1 = (cupy.cumsum(diff_y1 * input_slice) +
starting_diff[0] * input[1] +
starting_diff[1] * input[0])
y1 = _find_initial_cond(
all_valid[1:][:input_slice.size], cum_poly_y1, input.size, i)
if not cupy.any(cupy.isnan(cupy.r_[y0, y1])):
break
if cupy.any(cupy.isnan(cupy.r_[y0, y1])):
raise ValueError(
'Sum to find symmetric boundary conditions did not converge.')
# Apply the system cs / (1 - a2 * z^-1 - a3 * z^-2)
zi = cupy.r_[y0, y1]
y_fwd, _ = lfilter(cs, cupy.r_[1, -a2, -a3], input[2:], zi=zi)
y_fwd = cupy.r_[zi, y_fwd]
# Then compute the symmetric backward starting conditions
compute_symiirorder2_bwd_sc = _get_module_func(
SYMIIR2_MODULE, 'compute_symiirorder2_bwd_sc', cs)
diff = cupy.empty((block_sz,), dtype=cs.dtype)
all_valid = cupy.empty((block_sz,), dtype=cupy.bool_)
rev_input = input[::-1]
y0 = cupy.nan
for i in range(0, input.size + 1, block_sz):
compute_symiirorder2_bwd_sc(
(1,), (block_sz,), (
input.size + 1, i, 0, 1, cs, cupy.asarray(rsq, cs.dtype),
cupy.asarray(omega, cs.dtype), precision, all_valid, diff))
input_slice = rev_input[i:i + block_sz]
cum_poly_y0 = cupy.cumsum(diff[:input_slice.size] * input_slice)
y0 = _find_initial_cond(
all_valid[:input_slice.size], cum_poly_y0, input.size, i)
if not cupy.isnan(y0):
break
if cupy.isnan(y0):
raise ValueError(
'Sum to find symmetric boundary conditions did not converge.')
y1 = cupy.nan
for i in range(0, input.size + 1, block_sz):
compute_symiirorder2_bwd_sc(
(1,), (block_sz,), (
input.size + 1, i, -1, 2, cs, cupy.asarray(rsq, cs.dtype),
cupy.asarray(omega, cs.dtype), precision, all_valid, diff))
input_slice = rev_input[i:i + block_sz]
cum_poly_y1 = cupy.cumsum(diff[:input_slice.size] * input_slice)
y1 = _find_initial_cond(
all_valid[:input_slice.size], cum_poly_y1, input.size, i)
if not cupy.isnan(y1):
break
if cupy.isnan(y1):
raise ValueError(
'Sum to find symmetric boundary conditions did not converge.')
# Apply the system cs / (1 - a2 * z^1 - a3 * z^2)
zi = cupy.r_[y0, y1]
out, _ = lfilter(cs, cupy.r_[1, -a2, -a3], y_fwd[:-2][::-1], zi=zi)
return cupy.r_[out[::-1], zi[::-1]]