-
Notifications
You must be signed in to change notification settings - Fork 7
/
Copy pathsupdeq_deq.m
executable file
·423 lines (362 loc) · 19 KB
/
supdeq_deq.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
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
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
%% SUpDEq - Spatial Upsampling by Directional Equalization
%
% function [denseHRTFdataset, denseHRIRdataset, denseHRTFdataset_sfd] = supdeq_deq(eqHRTFdataset, deqDataset, N, samplingGrid, windowLength, normalization, parallax)
%
% This function performs the de-equalization of a equalized sparse HRTF
% dataset, here called eqHRTFdataset, (in SH-domain / stored as SH-coefficients)
% with the pre-defined deqDataset (dataset for de-equalization in SH-domain /
% stored as SH-coefficients). The result is a dense hrtf dataset (according
% to the (dense) spatial sampling grid), saved as in structs as HRTFs
% (denseHRTFdataset), HRIRs (denseHRIRdataset) or as SH-coefficients
% (denseHRTFdataset_sfd)
%
% Output:
% denseHRTFdataset - Dense HRTF dataset as a struct
% denseHRIRdataset - Dense HRIR dataset as a struct
% denseHRTFdataset_sfd - Dense HRIR dataset as SH-coefficients,
% transformed with given order N
%
% Input:
% eqHRTFdataset - Struct with (SH-Coefficients) of the equalized
% HRTF dataset for the left (Hl_nm) and right (Hr_nm)
% channel/ear, absolute frequency scale f,
% transform order N, and FFToversize
% deqDataset - Struct with de-equalization dataset in SH-domain / as
% SH-coefficients. Can be the output of supdeq_getEqDataset.
% N - Transform order N
% samplingGrid - Spatial sampling grid (Q x 2 or Q x 3 matrix),
% where the first column holds the azimuth,
% the second the elevation (both in degree), and
% optionally the third the sampling weights.
% Azimuth in degree
% (0=front, 90=left, 180=back, 270=right)
% (0 points to positive x-axis, 90 to positive y-axis)
% Elevations in degree
% (0=North Pole, 90=front, 180=South Pole)
% (0 points to positive z-axis, 180 to negative z-axis)
% windowLength - Values of head and tail window in samples [a,b].
% Head and tail window will be a half sided hann
% window with the defined length. Windowing will be
% applied to the HRIRs. The HRTF dataset in frequency-
% and SH-domain will be calculated based on these
% windowed HRIRs. If the matrix is empty ([]), no
% windowing will be applied!
% Default: [], no windowing
% normalization - Linear value (for example 0.99) which can be
% applied in order to normalize the de-equalized
% HRIRs. Can be useful if de-equalization with
% waveType == 1 (spherical wave) is applied.
% Default: [], no normalization.
% parallax - Boolean with true or false. Only applies if
% waveType in deq is 1 (spherical wave) and a
% distance shift is applied (DVF).
% false: acoustic parallax not considered
% true (Default: acoustic parallax considered
%
% Dependencies: SOFiA toolbox
%
% (C) 2018 by JMA, Johannes M. Arend
% CP, Christoph Pörschmann
% TH Köln - University of Applied Sciences
% Institute of Communications Engineering
% Department of Acoustics and Audio Signal Processing
function [denseHRTFdataset, denseHRIRdataset, denseHRTFdataset_sfd] = supdeq_deq(eqHRTFdataset, deqDataset, N, samplingGrid, windowLength, normalization, parallax)
if nargin < 5 || isempty(windowLength)
windowLength = [];
end
if nargin <6 || isempty(normalization)
normalization = [];
end
if nargin < 7 || isempty(parallax)
parallax = true;
end
%Check if eqDataset was with limiting but deqDataset is without limiting
if isfield(eqHRTFdataset,'limited')
if eqHRTFdataset.limited
if isfield(deqDataset,'limited')
if ~deqDataset.limited
error('Limited eq transfer functions applied for equalization, but not applied for de-equalization!');
end
else
error('Limited eq transfer functions applied for equalization, but not applied for de-equalization!');
end
end
end
%Check if applied radius and fA is the same
if isfield(eqHRTFdataset,'limited') && isfield(deqDataset,'limited')
if eqHRTFdataset.limited && deqDataset.limited
if eqHRTFdataset.limitInfo.appliedRadius ~= deqDataset.limitInfo.appliedRadius
error('Applied radius of limited equalization and de-equalization dataset must be the same!');
end
if eqHRTFdataset.limitInfo.fA ~= deqDataset.limitInfo.fA
error('fA of limited equalization and de-equalization dataset must be the same!');
end
end
end
%Check if phaseOnly equalization was applied
if isfield(eqHRTFdataset,'phaseOnly')
phaseOnly = 1;
else
phaseOnly = 0;
end
%Get fs
fs = eqHRTFdataset.f(end)*2;
%% PROCESSING
%Processing for de-equalization with plane wave deqDataset
if deqDataset.waveType == 0
%Transform deqDataset to Fourier domain at dense sampling grid points
%(inverse spherical Fourier transform)
fprintf('Extracting %d deq transfer functions. This may take some time...\n',size(samplingGrid,1));
[deqTF_L,deqTF_R] = supdeq_getArbHRTF(deqDataset,samplingGrid,[],[],'ak'); %Use AKtools...
fprintf('%d deq transfer functions extracted...\n',size(samplingGrid,1))
if phaseOnly
fprintf('Phase only de-equalization...\n')
%Get only phase response of eqTF
deqTF_L_phase = angle(deqTF_L);
deqTF_R_phase = angle(deqTF_R);
deqTF_L_mag = ones(size(deqTF_L,1),size(deqTF_L,2));
deqTF_R_mag = ones(size(deqTF_R,1),size(deqTF_R,2));
deqTF_L = deqTF_L_mag.*exp(1i*deqTF_L_phase);
deqTF_R = deqTF_R_mag.*exp(1i*deqTF_R_phase);
end
%Get HRTF from equalized HRTF dataset by inverse spherical Fourier transform at
%sampling points of dense sampling grid
fprintf('Extracting %d HRTFs. This may take some time...\n',size(samplingGrid,1));
[eqHRTF_L,eqHRTF_R] = supdeq_getArbHRTF(eqHRTFdataset,samplingGrid,[],[],'ak'); %Use AKtools...
fprintf('%d HRTFs extracted...\n',size(samplingGrid,1))
%Perform de-equalization (multiplication in frequency domain)
HRTF_deequalized_L = zeros(size(samplingGrid,1),length(deqDataset.f));
HRTF_deequalized_R = zeros(size(samplingGrid,1),length(deqDataset.f));
for kk = 1:size(samplingGrid,1)
HRTF_deequalized_L(kk,:)=eqHRTF_L(kk,:).*deqTF_L(kk,:);
HRTF_deequalized_R(kk,:)=eqHRTF_R(kk,:).*deqTF_R(kk,:);
if ~mod(kk,100)
fprintf('%d of %d HRTFs de-equalized...\n',kk,size(samplingGrid,1))
end
end
end
%Processing for de-equalization with spherical wave deqDataset (DVF - Distance Variation Function)
if deqDataset.waveType == 1
fprintf('Applying DVF...\n');
%Transform deqDataset to Fourier domain at dense sampling grid points
%(inverse spherical Fourier transform)
fprintf('Extracting %d deq transfer functions. This may take some time...\n',size(samplingGrid,1));
[deqTF_L,deqTF_R] = supdeq_getArbHRTF(deqDataset,samplingGrid,[],[],'ak'); %Use AKtools...
fprintf('%d deq transfer functions extracted...\n',size(samplingGrid,1))
if parallax %Consider acoustic parallax
fprintf('Considering acoustic parallax...\n');
%Get parallax shifted dense sampling grid for left and right ear
[samplingGridParL,samplingGridParR] = supdeq_parallax(samplingGrid,deqDataset.sourceDistance,deqDataset.radius,deqDataset.earPosition);
%Get HRTF from equalized HRTF dataset by inverse spherical Fourier transform at
%sampling points of dense sampling grid
fprintf('Extracting %d HRTFs. This may take some time...\n',size(samplingGridParL,1));
%Get left and right HRTFs separately
[eqHRTF_L,~] = supdeq_getArbHRTF(eqHRTFdataset,samplingGridParL,'DEG',0,'ak'); %Use AKtools...
[~,eqHRTF_R] = supdeq_getArbHRTF(eqHRTFdataset,samplingGridParR,'DEG',1,'ak'); %Use AKtools...
fprintf('%d HRTFs extracted...\n',size(samplingGrid,1))
else %parallax = false, do not consider acoustic parallax
fprintf('Not considering acoustic parallax...\n');
%Get HRTF from equalized HRTF dataset by inverse spherical Fourier transform at
%sampling points of dense sampling grid
fprintf('Extracting %d HRTFs. This may take some time...\n',size(samplingGrid,1));
[eqHRTF_L,eqHRTF_R] = supdeq_getArbHRTF(eqHRTFdataset,samplingGrid,[],[],'ak'); %Use AKtools...
fprintf('%d HRTFs extracted...\n',size(samplingGrid,1))
end
%Perform de-equalization (multiplication in frequency domain)
HRTF_deequalized_L = zeros(size(samplingGrid,1),length(deqDataset.f));
HRTF_deequalized_R = zeros(size(samplingGrid,1),length(deqDataset.f));
for kk = 1:size(samplingGrid,1)
HRTF_deequalized_L(kk,:)=eqHRTF_L(kk,:).*deqTF_L(kk,:);
HRTF_deequalized_R(kk,:)=eqHRTF_R(kk,:).*deqTF_R(kk,:);
if ~mod(kk,100)
fprintf('%d of %d HRTFs de-equalized...\n',kk,size(samplingGrid,1))
end
end
end
%% Store in structs
% Transform HRTFs to time domain --> get HRIRs
%Get both sided spectrum
HRTF_deequalized_L_double = conj([HRTF_deequalized_L(:,:), conj(fliplr(HRTF_deequalized_L(:,2:end-1)))])';
HRTF_deequalized_R_double = conj([HRTF_deequalized_R(:,:), conj(fliplr(HRTF_deequalized_R(:,2:end-1)))])';
%Transform back to time domain
HRIR_deequalized_L = real(ifft(HRTF_deequalized_L_double));
HRIR_deequalized_R = real(ifft(HRTF_deequalized_R_double));
%Cut if FFToversize of eqHRTFdataset > 1
if eqHRTFdataset.FFToversize > 1
newLength = size(HRIR_deequalized_L,1)/eqHRTFdataset.FFToversize;
HRIR_deequalized_L = HRIR_deequalized_L(1:newLength,:);
HRIR_deequalized_R = HRIR_deequalized_R(1:newLength,:);
end
%Window if windowLength is not empty
if ~isempty(windowLength) %windowing
disp('Applying windows...');
%Apply normalization if desired
if ~isempty(normalization)
disp('Applying normalization...');
%Get abs max of left/right channel
maxL = max(max(abs(HRIR_deequalized_L)));
maxR = max(max(abs(HRIR_deequalized_R)));
%Calculate normalization factor
if maxL >= maxR
normFactor = normalization/maxL;
elseif maxR > maxL
normFactor = normalization/maxR;
end
%Apply to HRIRs
HRIR_deequalized_L = HRIR_deequalized_L*normFactor;
HRIR_deequalized_R = HRIR_deequalized_R*normFactor;
end
%Windowing
%Get window function in time domain
win = supdeq_win(size(HRIR_deequalized_L,1),windowLength);
%Apply to HRIRs
for kk = 1:size(HRIR_deequalized_L,2)
HRIR_deequalized_L(:,kk) = HRIR_deequalized_L(:,kk) .* win;
HRIR_deequalized_R(:,kk) = HRIR_deequalized_R(:,kk) .* win;
end
%Store in struct denseHRIRdataset (window, normalization optionally)
denseHRIRdataset.HRIR_L = HRIR_deequalized_L'; %Flip again to have same array settings as in SOFiA
denseHRIRdataset.HRIR_R = HRIR_deequalized_R';
denseHRIRdataset.fs = eqHRTFdataset.f(end)*2;
denseHRIRdataset.samplingGrid = samplingGrid;
denseHRIRdataset.deqEarDistance = deqDataset.earDistance;
denseHRIRdataset.deqWaveType = deqDataset.waveType;
if deqDataset.waveType == 1
denseHRIRdataset.sourceDistance = deqDataset.sourceDistance;
end
if ~isempty(normalization)
denseHRIRdataset.normFactor = normFactor;
end
denseHRIRdataset.windowLength = windowLength;
%Transform de-equalized windowed denseHRIRdataset to Fourier domain to get
%denseHRTFdataset (with same FFToversize)
HRTF_deequalized_win_L = fft(HRIR_deequalized_L',size(HRIR_deequalized_L,1)*eqHRTFdataset.FFToversize,2);
HRTF_deequalized_win_L = HRTF_deequalized_win_L(:,1:end/2+1);
HRTF_deequalized_win_R = fft(HRIR_deequalized_R',size(HRIR_deequalized_R,1)*eqHRTFdataset.FFToversize,2);
HRTF_deequalized_win_R = HRTF_deequalized_win_R(:,1:end/2+1);
%Store in struct denseHRTFdataset (window, normalization optionally)
denseHRTFdataset.HRTF_L = HRTF_deequalized_win_L;
denseHRTFdataset.HRTF_R = HRTF_deequalized_win_R;
denseHRTFdataset.f = eqHRTFdataset.f;
denseHRTFdataset.fs = eqHRTFdataset.f(end)*2;
denseHRTFdataset.FFToversize = eqHRTFdataset.FFToversize;
denseHRTFdataset.samplingGrid = samplingGrid;
denseHRTFdataset.deqEarDistance = deqDataset.earDistance;
denseHRTFdataset.deqWaveType = deqDataset.waveType;
if deqDataset.waveType == 1
denseHRTFdataset.sourceDistance = deqDataset.sourceDistance;
end
if ~isempty(normalization)
denseHRTFdataset.normFactor = normFactor;
end
denseHRTFdataset.windowLength = windowLength;
% Transform de-equalized HRTFs to SH-domain --> get
% denseHRTFdataset_sfd (window, normalization optionally)
fprintf('Performing spherical Fourier transform with N = %d. This may take some time...\n',N);
denseHRTFdataset_sfd = supdeq_hrtf2sfd(HRTF_deequalized_win_L,HRTF_deequalized_win_R,N,samplingGrid,fs,'ak'); %Use AKtools...
denseHRTFdataset_sfd.FFToversize = eqHRTFdataset.FFToversize;
denseHRTFdataset_sfd.deqEarDistance = deqDataset.earDistance;
denseHRTFdataset_sfd.deqWaveType = deqDataset.waveType;
if deqDataset.waveType == 1
denseHRTFdataset_sfd.sourceDistance = deqDataset.sourceDistance;
end
if ~isempty(normalization)
denseHRTFdataset_sfd.normFactor = normFactor;
end
denseHRTFdataset_sfd.windowLength = windowLength;
fprintf('Done with de-equalization...\n');
else %no windowing
%Apply normalization if desired
if ~isempty(normalization)
disp('Applying normalization...');
%Get abs max of left/right channel
maxL = max(max(abs(HRIR_deequalized_L)));
maxR = max(max(abs(HRIR_deequalized_R)));
%Calculate normalization factor
if maxL >= maxR
normFactor = normalization/maxL;
elseif maxR > maxL
normFactor = normalization/maxR;
end
%Apply to HRIRs
HRIR_deequalized_L = HRIR_deequalized_L*normFactor;
HRIR_deequalized_R = HRIR_deequalized_R*normFactor;
%Store in struct denseHRIRdataset (no window, normalization)
denseHRIRdataset.HRIR_L = HRIR_deequalized_L'; %Flip again to have same array settings as in SOFiA
denseHRIRdataset.HRIR_R = HRIR_deequalized_R';
denseHRIRdataset.fs = eqHRTFdataset.f(end)*2;
denseHRIRdataset.samplingGrid = samplingGrid;
denseHRIRdataset.deqEarDistance = deqDataset.earDistance;
denseHRIRdataset.deqWaveType = deqDataset.waveType;
if deqDataset.waveType == 1
denseHRIRdataset.sourceDistance = deqDataset.sourceDistance;
end
denseHRIRdataset.normFactor = normFactor;
%Transform de-equalized normalized denseHRIRdataset to Fourier domain to get
%denseHRTFdataset (with same FFToversize)
HRTF_deequalized_L = fft(HRIR_deequalized_L',size(HRIR_deequalized_L,1)*eqHRTFdataset.FFToversize,2);
HRTF_deequalized_L = HRTF_deequalized_L(:,1:end/2+1);
HRTF_deequalized_R = fft(HRIR_deequalized_R',size(HRIR_deequalized_R,1)*eqHRTFdataset.FFToversize,2);
HRTF_deequalized_R = HRTF_deequalized_R(:,1:end/2+1);
%Store in struct denseHRTFdataset (no window, normalization)
denseHRTFdataset.HRTF_L = HRTF_deequalized_L;
denseHRTFdataset.HRTF_R = HRTF_deequalized_R;
denseHRTFdataset.f = eqHRTFdataset.f;
denseHRTFdataset.fs = eqHRTFdataset.f(end)*2;
denseHRTFdataset.FFToversize = eqHRTFdataset.FFToversize;
denseHRTFdataset.samplingGrid = samplingGrid;
denseHRTFdataset.deqEarDistance = deqDataset.earDistance;
denseHRTFdataset.deqWaveType = deqDataset.waveType;
if deqDataset.waveType == 1
denseHRTFdataset.sourceDistance = deqDataset.sourceDistance;
end
denseHRTFdataset.normFactor = normFactor;
% Transform de-equalized HRTFs to SH-domain --> get
% denseHRTFdataset_sfd (no window, normalization)
fprintf('Performing spherical Fourier transform with N = %d. This may take some time...\n',N);
denseHRTFdataset_sfd = supdeq_hrtf2sfd(HRTF_deequalized_L,HRTF_deequalized_R,N,samplingGrid,fs,'ak'); %Use AKtools...
denseHRTFdataset_sfd.FFToversize = eqHRTFdataset.FFToversize;
denseHRTFdataset_sfd.deqEarDistance = deqDataset.earDistance;
denseHRTFdataset_sfd.deqWaveType = deqDataset.waveType;
if deqDataset.waveType == 1
denseHRTFdataset_sfd.sourceDistance = deqDataset.sourceDistance;
end
denseHRTFdataset_sfd.normFactor = normFactor;
fprintf('Done with de-equalization...\n');
else % No normalization
%Store in struct denseHRIRdataset (no window, no normalization)
denseHRIRdataset.HRIR_L = HRIR_deequalized_L'; %Flip again to have same array settings as in SOFiA
denseHRIRdataset.HRIR_R = HRIR_deequalized_R';
denseHRIRdataset.fs = eqHRTFdataset.f(end)*2;
denseHRIRdataset.samplingGrid = samplingGrid;
denseHRIRdataset.deqEarDistance = deqDataset.earDistance;
denseHRIRdataset.deqWaveType = deqDataset.waveType;
if deqDataset.waveType == 1
denseHRIRdataset.sourceDistance = deqDataset.sourceDistance;
end
%Store in struct denseHRTFdataset (no window, no normalization)
denseHRTFdataset.HRTF_L = HRTF_deequalized_L;
denseHRTFdataset.HRTF_R = HRTF_deequalized_R;
denseHRTFdataset.f = eqHRTFdataset.f;
denseHRTFdataset.fs = eqHRTFdataset.f(end)*2;
denseHRTFdataset.FFToversize = eqHRTFdataset.FFToversize;
denseHRTFdataset.samplingGrid = samplingGrid;
denseHRTFdataset.deqEarDistance = deqDataset.earDistance;
denseHRTFdataset.deqWaveType = deqDataset.waveType;
if deqDataset.waveType == 1
denseHRTFdataset.sourceDistance = deqDataset.sourceDistance;
end
% Transform de-equalized HRTFs to SH-domain --> get
% denseHRTFdataset_sfd (no window, no normalization)
fprintf('Performing spherical Fourier transform with N = %d. This may take some time...\n',N);
denseHRTFdataset_sfd = supdeq_hrtf2sfd(HRTF_deequalized_L,HRTF_deequalized_R,N,samplingGrid,fs,'ak'); %Use AKtools...
denseHRTFdataset_sfd.FFToversize = eqHRTFdataset.FFToversize;
denseHRTFdataset_sfd.deqEarDistance = deqDataset.earDistance;
denseHRTFdataset_sfd.deqWaveType = deqDataset.waveType;
if deqDataset.waveType == 1
denseHRTFdataset_sfd.sourceDistance = deqDataset.sourceDistance;
end
fprintf('Done with de-equalization...\n');
end
end
end