forked from sfstoolbox/sfs-matlab
-
Notifications
You must be signed in to change notification settings - Fork 0
/
pwd_imp_circexp.m
72 lines (66 loc) · 3.31 KB
/
pwd_imp_circexp.m
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
function ppwd = pwd_imp_circexp(pm,Npw)
%PWD_IMP_CIRCEXP converts a circular basis expansion of a sound field to its
%two-dimensional plane wave decomposition
%
% Usage: ppwd = pwd_imp_circexp(pm,[Npw])
%
% Input parameters:
% pm - circular basis expansion [N x (M+1)]
% Npw - number of equi-angular distributed plane waves, optional,
% default: 2*M+1
%
% Output parameters:
% ppwd - plane wave decomposition [N x Npw]
%
% See also: driving_function_imp_localwfs_sbl_ps,
% driving_function_imp_localwfs_sbl_pw
%*****************************************************************************
% The MIT License (MIT) *
% *
% Copyright (c) 2010-2019 SFS Toolbox Developers *
% *
% 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. *
% *
% The SFS Toolbox allows to simulate and investigate sound field synthesis *
% methods like wave field synthesis or higher order ambisonics. *
% *
% https://sfs.readthedocs.io [email protected] *
%*****************************************************************************
%% ===== Checking of input parameters ==================================
nargmin = 1;
nargmax = 2;
narginchk(nargmin,nargmax);
isargmatrix(pm);
M = size(pm,2)-1;
if nargin == nargmin
Npw = 2*M+1;
else
isargpositivescalar(Npw);
end
%% ===== Computation ====================================================
% Implementation of
% ___
% _ \
% p(phipw, t) = /__ p (t) i^m e^(-i m phipw)
% m=-M..M m
% with
%
% phipw = n * 2*pi/Npw
pm = [conj(pm(:,end:-1:2)), pm]; % append coefficients for negative m
ppwd = inverse_cht(bsxfun(@times,pm,1i.^(-M:M)),Npw);