/
_convolve.py
223 lines (185 loc) · 7.03 KB
/
_convolve.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
"""
Some of the functions defined here were ported directly from CuSignal under
terms of the MIT license, under the following notice:
Copyright (c) 2019-2023 NVIDIA CORPORATION & AFFILIATES. All rights reserved.
Permission is hereby granted, free of charge, to any person obtaining a
copy of this software and associated documentation files (the "Software"),
to deal in the Software without restriction, including without limitation
the rights to use, copy, modify, merge, publish, distribute, sublicense,
and/or sell copies of the Software, and to permit persons to whom the
Software is furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in
all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
DEALINGS IN THE SOFTWARE.
"""
import cupy
_convolve1d2o_kernel = cupy.ElementwiseKernel(
'raw T in1, raw T in2, int32 W, int32 H', 'T out',
"""
T temp {};
for (int x = 0; x < W; x++) {
for (int y = 0; y < H; y++) {
temp += in1[i + W - x - 1] * in1[i + H - y - 1] * in2[H * x + y];
}
}
out = temp;
""",
"cupy_convolved2o",
)
def _convolve1d2o(in1, in2, mode):
assert mode == "valid"
out_dim = in1.shape[0] - max(in2.shape) + 1
dtype = cupy.result_type(in1, in2)
out = cupy.empty(out_dim, dtype=dtype)
_convolve1d2o_kernel(in1, in2, *in2.shape, out)
return out
def convolve1d2o(in1, in2, mode='valid', method='direct'):
"""
Convolve a 1-dimensional arrays with a 2nd order filter.
This results in a second order convolution.
Convolve `in1` and `in2`, with the output size determined by the
`mode` argument.
Parameters
----------
in1 : array_like
First input.
in2 : array_like
Second input. Should have the same number of dimensions as `in1`.
mode : str {'full', 'valid', 'same'}, optional
A string indicating the size of the output:
``full``
The output is the full discrete linear convolution
of the inputs. (Default)
``valid``
The output consists only of those elements that do not
rely on the zero-padding. In 'valid' mode, either `in1` or `in2`
must be at least as large as the other in every dimension.
``same``
The output is the same size as `in1`, centered
with respect to the 'full' output.
method : str {'auto', 'direct', 'fft'}, optional
A string indicating which method to use to calculate the convolution.
``direct``
The convolution is determined directly from sums, the definition of
convolution.
``fft``
The Fourier Transform is used to perform the convolution by calling
`fftconvolve`.
``auto``
Automatically chooses direct or Fourier method based on an estimate
of which is faster (default).
Returns
-------
out : ndarray
A 1-dimensional array containing a subset of the discrete linear
convolution of `in1` with `in2`.
See Also
--------
convolve
convolve1d2o
convolve1d3o
Examples
--------
Convolution of a 2nd order filter on a 1d signal
>>> import cusignal as cs
>>> import numpy as np
>>> d = 50
>>> a = np.random.uniform(-1,1,(200))
>>> b = np.random.uniform(-1,1,(d,d))
>>> c = cs.convolve1d2o(a,b)
"""
if in1.ndim != 1:
raise ValueError('in1 should have one dimension')
if in2.ndim != 2:
raise ValueError('in2 should have three dimension')
if mode in ["same", "full"]:
raise NotImplementedError("Mode == {} not implemented".format(mode))
if method == "direct":
return _convolve1d2o(in1, in2, mode)
else:
raise NotImplementedError("Only Direct method implemented")
_convolve1d3o_kernel = cupy.ElementwiseKernel(
'raw T in1, raw T in2, int32 W, int32 H, int32 D', 'T out',
"""
T temp {};
for (int x = 0; x < W; x++) {
for (int y = 0; y < H; y++) {
for (int z = 0; z < D; z++) {
temp += in1[i + W - x - 1] * in1[i + H - y - 1] *
in1[i + D - z - 1] * in2[(H * x + y) * D + z];
}
}
}
out = temp;
""",
"cupy_convolved3o",
)
def _convolve1d3o(in1, in2, mode):
assert mode == "valid"
out_dim = in1.shape[0] - max(in2.shape) + 1
dtype = cupy.result_type(in1, in2)
out = cupy.empty(out_dim, dtype=dtype)
_convolve1d3o_kernel(in1, in2, *in2.shape, out)
return out
def convolve1d3o(in1, in2, mode='valid', method='direct'):
"""
Convolve a 1-dimensional array with a 3rd order filter.
This results in a third order convolution.
Convolve `in1` and `in2`, with the output size determined by the
`mode` argument.
Parameters
----------
in1 : array_like
First input. Should have one dimension.
in2 : array_like
Second input. Should have three dimensions.
mode : str {'full', 'valid', 'same'}, optional
A string indicating the size of the output:
``full``
The output is the full discrete linear convolution
of the inputs. (Default)
``valid``
The output consists only of those elements that do not
rely on the zero-padding. In 'valid' mode, either `in1` or `in2`
must be at least as large as the other in every dimension.
``same``
The output is the same size as `in1`, centered
with respect to the 'full' output.
method : str {'auto', 'direct', 'fft'}, optional
A string indicating which method to use to calculate the convolution.
``direct``
The convolution is determined directly from sums, the definition of
convolution.
``fft``
The Fourier Transform is used to perform the convolution by calling
`fftconvolve`.
``auto``
Automatically chooses direct or Fourier method based on an estimate
of which is faster (default).
Returns
-------
out : ndarray
A 1-dimensional array containing a subset of the discrete linear
convolution of `in1` with `in2`.
See Also
--------
convolve
convolve1d2o
convolve1d3o
"""
if in1.ndim != 1:
raise ValueError('in1 should have one dimension')
if in2.ndim != 3:
raise ValueError('in2 should have three dimension')
if mode in ["same", "full"]:
raise NotImplementedError("Mode == {} not implemented".format(mode))
if method == "direct":
return _convolve1d3o(in1, in2, mode)
else:
raise NotImplementedError("Only Direct method implemented")