Skip to content

Commit 9cdbcfa

Browse files
author
Till Rettberg
committed
Merge in circular radial functions
2 parents 30f6d3a + 1cf5efd commit 9cdbcfa

File tree

5 files changed

+347
-0
lines changed

5 files changed

+347
-0
lines changed
Lines changed: 48 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,48 @@
1+
"""
2+
Compute the plane wave decomposition for an incident broadband plane wave
3+
on an open circular array using a modal beamformer of finite order.
4+
"""
5+
6+
import numpy as np
7+
import matplotlib.pyplot as plt
8+
import micarray
9+
10+
N = 90 # order of modal beamformer/microphone array
11+
pw_angle = 1.23 * np.pi # incidence angle of plane wave
12+
pol_pwd = np.linspace(0, 2*np.pi, 91, endpoint=False) # angles for plane wave decomposition
13+
k = np.linspace(0.1, 20, 100) # wavenumber vector
14+
r = 1 # radius of array
15+
16+
# get uniform grid (microphone positions) of order N
17+
pol, weights = micarray.modal.angular.grid_equal_polar_angle(N)
18+
# get circular harmonics matrix for sensors
19+
Psi_p = micarray.modal.angular.cht_matrix(N, pol, weights)
20+
# get circular harmonics matrix for a source ensemble of azimuthal plane wave
21+
Psi_q = micarray.modal.angular.cht_matrix(N, pol_pwd)
22+
# get radial filters
23+
Bn = micarray.modal.radial.circular_pw(N, k, r, setup='open')
24+
Dn, _ = micarray.modal.radial.regularize(1/Bn, 100, 'softclip')
25+
D = micarray.modal.radial.circ_diagonal_mode_mat(Dn)
26+
27+
# compute microphone signals for an incident broad-band plane wave
28+
p = np.exp(1j * k[:, np.newaxis]*r * np.cos(pol - pw_angle))
29+
# compute plane wave decomposition
30+
A_pwd = np.matmul(np.matmul(np.conj(Psi_q.T), D), Psi_p)
31+
q_pwd = np.squeeze(np.matmul(A_pwd, np.expand_dims(p, 2)))
32+
q_pwd_t = np.fft.fftshift(np.fft.irfft(q_pwd, axis=0), axes=0)
33+
34+
# visualize plane wave decomposition (aka beampattern)
35+
plt.figure()
36+
plt.pcolormesh(k, pol_pwd/np.pi, micarray.util.db(q_pwd.T), vmin=-40)
37+
plt.colorbar()
38+
plt.xlabel(r'$kr$')
39+
plt.ylabel(r'$\phi / \pi$')
40+
plt.title('Plane wave docomposition by modal beamformer (frequency domain)')
41+
plt.savefig('modal_open_beamformer_pwd_fd.png')
42+
43+
plt.figure()
44+
plt.pcolormesh(range(2*len(k)-2), pol_pwd/np.pi, micarray.util.db(q_pwd_t.T), vmin=-40)
45+
plt.colorbar()
46+
plt.ylabel(r'$\phi / \pi$')
47+
plt.title('Plane wave docomposition by modal beamformer (time domain)')
48+
plt.savefig('modal_open_beamformer_pwd_td.png')
Lines changed: 52 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,52 @@
1+
"""
2+
Compute the plane wave decomposition for an incident broadband plane wave
3+
on an rigid circular array using a modal beamformer of finite order.
4+
"""
5+
6+
import numpy as np
7+
import matplotlib.pyplot as plt
8+
import micarray
9+
10+
Nsf = 50 # order of the incident sound field
11+
N = 30 # order of modal beamformer/microphone array
12+
pw_angle = 1 * np.pi # incidence angle of plane wave
13+
pol_pwd = np.linspace(0, 2*np.pi, 180, endpoint=False) # angles for plane wave decomposition
14+
k = np.linspace(0.1, 20, 100) # wavenumber vector
15+
r = 1 # radius of array
16+
17+
# get uniform grid (microphone positions) of order N
18+
pol, weights = micarray.modal.angular.grid_equal_polar_angle(N)
19+
20+
# pressure on the surface of a rigid cylinder for an incident plane wave
21+
bn = micarray.modal.radial.circular_pw(Nsf, k, r, setup='rigid')
22+
D = micarray.modal.radial.circ_diagonal_mode_mat(bn)
23+
Psi_p = micarray.modal.angular.cht_matrix(Nsf, pol, weights)
24+
Psi_pw = micarray.modal.angular.cht_matrix(Nsf, pw_angle)
25+
p = np.matmul(np.matmul(np.conj(Psi_pw.T), D), Psi_p)
26+
p = np.squeeze(p)
27+
28+
# plane wave decomposition using modal beamforming
29+
Psi_p = micarray.modal.angular.cht_matrix(N, pol)
30+
Psi_q = micarray.modal.angular.cht_matrix(N, pol_pwd)
31+
Bn = micarray.modal.radial.circular_pw(N, k, r, setup='rigid')
32+
Dn, _ = micarray.modal.radial.regularize(1/Bn, 100, 'softclip')
33+
D = micarray.modal.radial.circ_diagonal_mode_mat(Dn)
34+
A_pwd = np.matmul(np.matmul(np.conj(Psi_q.T), D), Psi_p)
35+
q_pwd = np.squeeze(np.matmul(A_pwd, np.expand_dims(p, 2)))
36+
q_pwd_t = np.fft.fftshift(np.fft.irfft(q_pwd, axis=0), axes=0)
37+
38+
# visualize plane wave decomposition (aka beampattern)
39+
plt.figure()
40+
plt.pcolormesh(k, pol_pwd/np.pi, micarray.util.db(q_pwd.T), vmin=-40)
41+
plt.colorbar()
42+
plt.xlabel(r'$kr$')
43+
plt.ylabel(r'$\phi / \pi$')
44+
plt.title('Plane wave docomposition by modal beamformer (frequency domain)')
45+
plt.savefig('modal_open_beamformer_pwd_fd.png')
46+
47+
plt.figure()
48+
plt.pcolormesh(range(2*len(k)-2), pol_pwd/np.pi, micarray.util.db(q_pwd_t.T), vmin=-40)
49+
plt.colorbar()
50+
plt.ylabel(r'$\phi / \pi$')
51+
plt.title('Plane wave docomposition by modal beamformer (time domain)')
52+
plt.savefig('modal_open_beamformer_pwd_td.png')

micarray/modal/angular.py

Lines changed: 63 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -91,6 +91,48 @@ def Legendre_matrix(N, ctheta):
9191
return Lmn
9292

9393

94+
def cht_matrix(N, pol, weights=None):
95+
r"""Matrix of circular harmonics up to order N for given angles.
96+
97+
Computes a matrix of circular harmonics up to order :math:`N`
98+
for the given angles/grid.
99+
100+
.. math::
101+
\Psi = \left[ \begin{array}{ccccccc}
102+
1 & e^{i\varphi[0]} & \cdots & e^{iN\varphi[0]} & e^{-iN\varphi[0]} & \cdots & e^{-i\varphi[0]} \\
103+
1 & e^{i\varphi[1]} & \cdots & e^{iN\varphi[1]} & e^{-iN\varphi[1]} & \cdots & e^{-i\varphi[1]} \\
104+
\vdots & \vdots & \vdots & \vdots & \vdots & \vdots & \vdots \\
105+
1 & e^{i\varphi[Q-1]} & \cdots & e^{iN\varphi[Q-1]} & e^{-iN\varphi[Q-1]} & \cdots & e^{-i\varphi[Q-1]}
106+
\end{array} \right]
107+
108+
Parameters
109+
----------
110+
N : int
111+
Maximum order.
112+
pol : (Q,) array_like
113+
Polar angle.
114+
weights : (Q,) array_like, optional
115+
Weights.
116+
117+
Returns
118+
-------
119+
Psi : (2N+1, Q) numpy.ndarray
120+
Matrix of circular harmonics.
121+
"""
122+
pol = util.asarray_1d(pol)
123+
if pol.ndim == 0:
124+
Q = 1
125+
else:
126+
Q = len(pol)
127+
if weights is None:
128+
weights = np.ones(Q)
129+
Psi = np.zeros([(2*N+1), Q], dtype=complex)
130+
order = np.roll(np.arange(-N, N+1), -N)
131+
for i, n in enumerate(order):
132+
Psi[i, :] = weights * np.exp(1j * n * pol)
133+
return Psi
134+
135+
94136
def grid_equal_angle(n):
95137
"""Equi-angular sampling points on a sphere.
96138
@@ -153,3 +195,24 @@ def grid_gauss(n):
153195
weights = np.repeat(weights, 2*n+2)
154196
weights *= np.pi / (n+1) # sum(weights) == 4pi
155197
return azi, elev, weights
198+
199+
200+
def grid_equal_polar_angle(n):
201+
"""Equi-angular sampling points on a circle.
202+
203+
Parameters
204+
----------
205+
n : int
206+
Maximum order
207+
208+
Returns
209+
-------
210+
pol : array_like
211+
Polar angle.
212+
weights : array_like
213+
Weights.
214+
"""
215+
num_mic = 2*n+1
216+
pol = np.linspace(0, 2*np.pi, num=num_mic, endpoint=False)
217+
weights = 1/num_mic * np.ones(num_mic)
218+
return pol, weights

micarray/modal/radial.py

Lines changed: 168 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -229,3 +229,171 @@ def repeat_n_m(v):
229229
"""
230230
krlist = [np.tile(v, (2*i+1, 1)).T for i, v in enumerate(v.T.tolist())]
231231
return np.squeeze(np.concatenate(krlist, axis=-1))
232+
233+
234+
def circular_pw(N, k, r, setup):
235+
r"""Radial coefficients for a plane wave.
236+
237+
Computes the radial component of the circular harmonics expansion of a
238+
plane wave impinging on a circular array.
239+
240+
.. math::
241+
242+
\mathring{P}_n(k) = i^{-n} b_n(kr)
243+
244+
Parameters
245+
----------
246+
N : int
247+
Maximum order.
248+
k : (M,) array_like
249+
Wavenumber.
250+
r : float
251+
Radius of microphone array.
252+
setup : {'open', 'card', 'rigid'}
253+
Array configuration (open, cardioids, rigid).
254+
255+
Returns
256+
-------
257+
bn : (M, N+1) numpy.ndarray
258+
Radial weights for all orders up to N and the given wavenumbers.
259+
"""
260+
kr = util.asarray_1d(k*r)
261+
n = np.arange(N+1)
262+
263+
bn = circ_radial_weights(N, kr, setup)
264+
for i, x in enumerate(kr):
265+
bn[i, :] = bn[i, :] * (1j)**(n)
266+
return bn
267+
268+
269+
def circular_ls(N, k, r, rs, setup):
270+
r"""Radial coefficients for a line source.
271+
272+
Computes the radial component of the circular harmonics expansion of a
273+
line source impinging on a circular array.
274+
275+
.. math::
276+
277+
\mathring{P}_n(k) = \frac{-i}{4} H_n^{(2)}(k r_s) b_n(kr)
278+
279+
Parameters
280+
----------
281+
N : int
282+
Maximum order.
283+
k : (M,) array_like
284+
Wavenumber.
285+
r : float
286+
Radius of microphone array.
287+
rs : float
288+
Distance of source.
289+
setup : {'open', 'card', 'rigid'}
290+
Array configuration (open, cardioids, rigid).
291+
292+
Returns
293+
-------
294+
bn : (M, N+1) numpy.ndarray
295+
Radial weights for all orders up to N and the given wavenumbers.
296+
"""
297+
k = util.asarray_1d(k)
298+
krs = k*rs
299+
n = np.arange(N+1)
300+
301+
bn = circ_radial_weights(N, k*r, setup)
302+
for i, x in enumerate(krs):
303+
Hn = special.hankel2(n, x)
304+
bn[i, :] = bn[i, :] * -1j/4 * Hn
305+
return bn
306+
307+
308+
def circ_radial_weights(N, kr, setup):
309+
r"""Radial weighing functions.
310+
311+
Computes the radial weighting functions for diferent array types
312+
313+
For instance for an rigid array
314+
315+
.. math::
316+
317+
b_n(kr) = J_n(kr) - \frac{J_n^\prime(kr)}{H_n^{(2)\prime}(kr)}H_n^{(2)}(kr)
318+
319+
Parameters
320+
----------
321+
N : int
322+
Maximum order.
323+
kr : (M,) array_like
324+
Wavenumber * radius.
325+
setup : {'open', 'card', 'rigid'}
326+
Array configuration (open, cardioids, rigid).
327+
328+
Returns
329+
-------
330+
bn : (M, N+1) numpy.ndarray
331+
Radial weights for all orders up to N and the given wavenumbers.
332+
333+
"""
334+
n = np.arange(N+1)
335+
Bns = np.zeros((len(kr), N+1), dtype=complex)
336+
for i, x in enumerate(kr):
337+
Jn = special.jv(n, x)
338+
if setup == 'open':
339+
bn = Jn
340+
elif setup == 'card':
341+
bn = Jn - 1j * special.jvp(n, x, n=1)
342+
elif setup == 'rigid':
343+
Jnd = special.jvp(n, x, n=1)
344+
Hn = special.hankel2(n, x)
345+
Hnd = special.h2vp(n, x)
346+
bn = Jn - Jnd/Hnd*Hn
347+
else:
348+
raise ValueError('setup must be either: open, card or rigid')
349+
Bns[i, :] = bn
350+
return np.squeeze(Bns)
351+
352+
353+
def circ_diagonal_mode_mat(bk):
354+
"""Diagonal matrix of radial coefficients for all modes/wavenumbers.
355+
356+
Parameters
357+
----------
358+
bk : (M, N+1) numpy.ndarray
359+
Vector containing values for all wavenumbers :math:`M` and modes up to
360+
order :math:`N`
361+
362+
Returns
363+
-------
364+
Bk : (M, 2*N+1, 2*N+1) numpy.ndarray
365+
Multidimensional array containing diagnonal matrices with input
366+
vector on main diagonal.
367+
"""
368+
bk = mirror_vec(bk)
369+
if len(bk.shape) == 1:
370+
bk = bk[np.newaxis, :]
371+
K, N = bk.shape
372+
Bk = np.zeros([K, N, N], dtype=complex)
373+
for k in range(K):
374+
Bk[k, :, :] = np.diag(bk[k, :])
375+
return np.squeeze(Bk)
376+
377+
378+
def mirror_vec(v):
379+
"""Mirror elements in a vector.
380+
381+
Returns a vector of length *2*len(v)-1* with symmetric elements.
382+
The first *len(v)* elements are the same as *v* and the last *len(v)-1*
383+
elements are *v[:0:-1]*. The function can be used to order the circular
384+
harmonic coefficients. If *v* is a matrix, it is treated as a stack of
385+
vectors residing in the last index and broadcast accordingly.
386+
387+
Parameters
388+
----------
389+
v : (, N+1) numpy.ndarray
390+
Input vector of stack of input vectors
391+
392+
Returns
393+
-------
394+
: (, 2*N+1) numpy.ndarray
395+
Vector of stack of vectors containing mirrored elements
396+
"""
397+
if len(v.shape) == 1:
398+
v = v[np.newaxis, :]
399+
return np.concatenate((v, v[:, :0:-1]), axis=1)

micarray/util.py

Lines changed: 16 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -91,3 +91,19 @@ def matdiagmul(A, b):
9191
for k in range(K):
9292
C[k, :, :] = A * b[k, :]
9393
return C
94+
95+
96+
def db(x, power=False):
97+
"""Convert *x* to decibel.
98+
99+
Parameters
100+
----------
101+
x : array_like
102+
Input data. Values of 0 lead to negative infinity.
103+
power : bool, optional
104+
If ``power=False`` (the default), *x* is squared before
105+
conversion.
106+
107+
"""
108+
with np.errstate(divide='ignore'):
109+
return 10 if power else 20 * np.log10(np.abs(x))

0 commit comments

Comments
 (0)