SAF
Loading...
Searching...
No Matches
saf_hades_synthesis.c
Go to the documentation of this file.
1/*
2 * This file is part of the saf_hades module.
3 * Copyright (c) 2021 - Leo McCormack & Janani Fernandez
4 *
5 * The saf_hades module is free software; you can redistribute it and/or
6 * modify it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or (at your
8 * option) any later version.
9 *
10 * The saf_hades module is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
12 * or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for
13 * more details.
14 *
15 * See <http://www.gnu.org/licenses/> for a copy of the GNU General Public
16 * License.
17 */
18
38#include "saf_hades_synthesis.h"
39#include "saf_hades_internal.h"
40
41#ifdef SAF_ENABLE_HADES_MODULE
42
43/* ========================================================================== */
44/* HADES Radial Editor */
45/* ========================================================================== */
46
62
64(
66)
67{
69
70 if (e != NULL) {
71 free(e);
72 e = NULL;
73 (*phREd) = NULL;
74 }
75}
76
78(
81 float dirGain_dB[360]
82)
83{
86 int band, edit_idx;
87 float azi, gain_lin;
88
89 for(band=0; band<e->nBands; band++){
90 /* Determine edit index */
91 azi = e->pGrid_dirs_deg[pcon->gains_idx[band]*2];
92 azi = azi < 0.0f ? azi+360.0f : azi; /* convert -180..180 if needed */
93 edit_idx = SAF_CLAMP((int)(azi + 0.5f), 0, 359); /* round to nearest integer */
94
95 /* Extra gain factor for the direct stream */
96 gain_lin = powf(10.0f, SAF_CLAMP(dirGain_dB[edit_idx], -60.0f, 12.0f)/20.0f);
97 pcon->gains_dir[band] *= gain_lin;
98 }
99}
100
101
102/* ========================================================================== */
103/* HADES Synthesis */
104/* ========================================================================== */
105
107(
108 hades_synthesis_handle* const phSyn,
109 hades_analysis_handle const hAna,
110 HADES_BEAMFORMER_TYPE beamOption,
111 int enableCM,
112 int refIndices[2],
113 hades_binaural_config* binConfig,
114 HADES_HRTF_INTERP_OPTIONS interpOption
115)
116{
118 *phSyn = (hades_synthesis_handle)s;
120 int band;
121 float_complex* H_W;
122 hades_binaural_config* bConfig;
123 const float_complex calpha = cmplxf(1.0f, 0.0f); const float_complex cbeta = cmplxf(0.0f, 0.0f); /* blas */
124
125 /* User configuration parameters */
126 s->beamOption = beamOption;
127 s->enableCM = enableCM;
128 s->refIndices[0] = refIndices[0];
129 s->refIndices[1] = refIndices[1];
130 s->interpOption = interpOption;
131
132 /* Default user parameters */
133 s->eq = malloc1d(a->nBands * sizeof(float));
134 s->streamBalance = malloc1d(a->nBands * sizeof(float));
135 for(band = 0; band<a->nBands; band++){
136 s->eq[band] = 1.0; /* Flat EQ */
137 s->streamBalance[band] = 1.0f; /* 50/50 direct/ambient balance (i.e., no biasing) */
138 }
139 s->synAvgCoeff = 1.0f - 1.0f/(4096.0f/a->blocksize); /* How much averaging of current mixing matrices with the previous mixing matrices */
140
141 /* Things relevant to the synthesiser, which are copied from the analyser to keep things aligned */
142 s->fbOpt = a->fbOpt;
143 s->nBands = a->nBands;
144 s->hopsize = a->hopsize;
145 s->blocksize = a->blocksize;
146 s->nGrid = a->nGrid;
147 s->nMics = a->nMics;
148 s->H_array = malloc1d(s->nBands*(s->nMics)*(s->nGrid)*sizeof(float_complex));
149 memcpy(s->H_array, a->H_array, s->nBands*(s->nMics)*(s->nGrid)*sizeof(float_complex));
150 s->DCM_array = malloc1d(s->nBands*(s->nMics)*(s->nMics)*sizeof(float_complex));
151 memcpy(s->DCM_array, a->DCM_array, s->nBands*(s->nMics)*(s->nMics)*sizeof(float_complex));
152 s->W = malloc1d(s->nGrid*(s->nGrid)*sizeof(float_complex));
153 memcpy(s->W, a->W, s->nGrid*(s->nGrid)*sizeof(float_complex));
154 s->grid_dirs_deg = malloc1d((s->nGrid)*2*sizeof(float));
155 memcpy(s->grid_dirs_deg, a->grid_dirs_deg, (s->nGrid)*2*sizeof(float));
156 s->grid_dirs_xyz = (float**)malloc2d((s->nGrid), 3, sizeof(float));
157 memcpy(FLATTEN2D(s->grid_dirs_xyz), a->grid_dirs_xyz, (s->nGrid)*3*sizeof(float));
158 s->timeSlots = a->timeSlots;
159 s->freqVector = malloc1d(s->nBands*sizeof(float));
160 memcpy(s->freqVector, a->freqVector, s->nBands*sizeof(float));
161
162 /* Time-frequency transform */
163 switch(s->fbOpt){
166 }
167
168 /* Copy binaural configuration */
170 bConfig = s->binConfig;
171 bConfig->lHRIR = binConfig->lHRIR;
172 bConfig->nHRIR = binConfig->nHRIR;
173 bConfig->hrir_fs = binConfig->hrir_fs;
174 bConfig->hrirs = malloc1d(bConfig->nHRIR * NUM_EARS * (bConfig->lHRIR) * sizeof(float));
175 memcpy(bConfig->hrirs, binConfig->hrirs, bConfig->nHRIR * NUM_EARS * (bConfig->lHRIR) * sizeof(float));
176 bConfig->hrir_dirs_deg = malloc1d(bConfig->nHRIR*2*sizeof(float));
177 memcpy(bConfig->hrir_dirs_deg, binConfig->hrir_dirs_deg, bConfig->nHRIR*2*sizeof(float));
178
179 /* Pre-process HRTFs, interpolate them for the scanning grid */
180 s->H_bin = calloc1d(s->nBands*NUM_EARS*(s->nGrid),sizeof(float_complex));
181 hades_getInterpolatedHRTFs(hAna, interpOption, bConfig, a->grid_dirs_deg, s->nGrid, s->H_bin);
182
183 /* Diffuse rendering variables */
184 s->DCM_bin_norm = malloc1d(a->nBands*NUM_EARS*NUM_EARS*sizeof(float_complex));
185 H_W = malloc1d(NUM_EARS*(a->nGrid)*sizeof(float_complex));
186 s->diffEQ = malloc1d(a->nBands*sizeof(float));
187 for(band=0; band<s->nBands; band++){
188 /* Binaural diffuse coherence matrix (not normalised yet!) */
189 cblas_cgemm(CblasRowMajor, CblasNoTrans, CblasNoTrans, NUM_EARS, a->nGrid, a->nGrid, &calpha,
190 &(s->H_bin[band*NUM_EARS*(a->nGrid)]), a->nGrid,
191 a->W, a->nGrid, &cbeta,
192 H_W, a->nGrid);
193 cblas_cgemm(CblasRowMajor, CblasNoTrans, CblasConjTrans, NUM_EARS, NUM_EARS, a->nGrid, &calpha,
194 H_W, a->nGrid,
195 &(s->H_bin[band*NUM_EARS*(a->nGrid)]), a->nGrid, &cbeta,
197 cblas_sscal(/*re+im*/2*NUM_EARS*NUM_EARS, 1.0f/(float)a->nGrid, (float*)&(s->DCM_bin_norm[band*NUM_EARS*NUM_EARS]), 1);
198
199 /* Compute EQ required to bring the overall diffuse-field magnitude response of the array to that of the HRTFs instead: */
200 /* sqrt(trace(H_bin_dcm(:,:,band))/(trace(H_grid_dcm(refIndices,refIndices,band))+eps)) */
201 s->diffEQ[band] = sqrtf((crealf(s->DCM_bin_norm[band*NUM_EARS*NUM_EARS]) + crealf(s->DCM_bin_norm[band*NUM_EARS*NUM_EARS + 3])) /
202 (crealf(s->DCM_array[band*(s->nMics)*(s->nMics) + s->refIndices[0]*(s->nMics) + s->refIndices[0]]) +
203 crealf(s->DCM_array[band*(s->nMics)*(s->nMics) + s->refIndices[1]*(s->nMics) + s->refIndices[1]]) + 2.23e-10f));
204 s->diffEQ[band] = SAF_MIN(s->diffEQ[band], 3.0f); /* Cap at a maximum of +9dB */
205
206 /* Normalise the binaural diffuse coherence matrix */
207 cblas_sscal(/*re+im*/2*NUM_EARS*NUM_EARS, 1.0f/(crealf(s->DCM_bin_norm[band*NUM_EARS*NUM_EARS]) + crealf(s->DCM_bin_norm[band*NUM_EARS*NUM_EARS + 3]) + 2.23e-10f),
208 (float*)&(s->DCM_bin_norm[band*NUM_EARS*NUM_EARS]), 1);
209 }
210 free(H_W);
211
212 /* Run-time variables */
213 utility_cpinv_create(&(s->hPinv), s->nMics, s->nMics);
216 s->As = malloc1d(s->nMics*sizeof(float_complex));
217 s->As_l = malloc1d(s->nMics*sizeof(float_complex));
218 s->As_r = malloc1d(s->nMics*sizeof(float_complex));
219 s->Q_diff = malloc1d(NUM_EARS*(s->nMics)*sizeof(float_complex));
220 s->Q_dir = malloc1d(NUM_EARS*(s->nMics)*sizeof(float_complex));
221 s->Q = malloc1d(NUM_EARS*(s->nMics)*sizeof(float_complex));
222 s->Cy = malloc1d(NUM_EARS*NUM_EARS*sizeof(float_complex));
223 s->new_M = malloc1d(NUM_EARS*(s->nMics)*sizeof(float_complex));
224 s->M = (float_complex**)malloc2d(s->nBands, NUM_EARS*(s->nMics), sizeof(float_complex));
225
226 /* Run-time audio buffers */
227 s->outTF = (float_complex***)malloc3d(s->nBands, NUM_EARS, s->timeSlots, sizeof(float_complex));
228 s->outTD = (float**)malloc2d(NUM_EARS, s->blocksize, sizeof(float));
229
230 /* Flush run-time buffers with zeros */
231 hades_synthesis_reset((*phSyn));
232}
233
235(
236 hades_synthesis_handle* const phSyn
237)
238{
240
241 if (s != NULL) {
242 /* Free user parameters */
243 free(s->eq);
244 free(s->streamBalance);
245 free(s->binConfig);
246
247 /* Free things copied from analyser */
248 free(s->H_array);
249 free(s->DCM_array);
250 free(s->W);
251 free(s->grid_dirs_deg);
252 free(s->grid_dirs_xyz);
253 free(s->freqVector);
254
255 /* Free time-frequency transform */
256 switch(s->fbOpt){
257 case HADES_USE_AFSTFT_LD: /* fall through */
258 case HADES_USE_AFSTFT: afSTFT_destroy(&(s->hFB_dec)); break;
259 }
260
261 /* HRTF and diffuse rendering variables */
262 free(s->H_bin);
263 free(s->DCM_bin_norm);
264 free(s->diffEQ);
265
266 /* Run-time variables */
270 free(s->As);
271 free(s->As_l);
272 free(s->As_r);
273 free(s->Q_diff);
274 free(s->Q_dir);
275 free(s->Q);
276 free(s->Cy);
277 free(s->new_M);
278 free(s->M);
279
280 /* Run-time audio buffers */
281 free(s->outTF);
282 free(s->outTD);
283
284 free(s);
285 s = NULL;
286 (*phSyn) = NULL;
287 }
288}
289
291(
292 hades_synthesis_handle const hSyn
293)
294{
296 if(hSyn==NULL)
297 return;
298 s = (hades_synthesis_data*)(hSyn);
299
300 /* Zero buffers, matrices etc. */
301 switch(s->fbOpt){
302 case HADES_USE_AFSTFT_LD: /* fall through */
304 }
305 memset(FLATTEN2D(s->M), 0, s->nBands*NUM_EARS*(s->nMics)*sizeof(float_complex));
306}
307
309(
310 hades_synthesis_handle const hSyn,
313 int nChannels,
314 int blocksize,
315 float** output
316)
317{
321 int i, j, ch, nMics, band, doa_idx, gain_idx;
322 float a, b, diffuseness, synAvgCoeff, streamBalance, eq, gain_dir, gain_diff, trace_M, reg_M, sum_As, targetEnergy;
323 float_complex g_l, g_r, h_dir[NUM_EARS], AsH_invCx_As;
324 float_complex Cx[HADES_MAX_NMICS*HADES_MAX_NMICS], conj_As[HADES_MAX_NMICS], AsH_invCx[HADES_MAX_NMICS*HADES_MAX_NMICS];
325 const float_complex calpha = cmplxf(1.0f, 0.0f); const float_complex cbeta = cmplxf(0.0f, 0.0f); /* blas */
326
327 nMics = s->nMics;
328 synAvgCoeff = SAF_CLAMP((s->synAvgCoeff), 0.0f, 0.99f);
329
330 /* Loop over bands and compute the mixing matrices */
331 for (band = 0; band < s->nBands; band++) {
332 /* Pull estimated (and possibly modified) spatial parameters for this band */
333 diffuseness = pcon->diffuseness[band];
334 saf_assert(diffuseness>-0.0001f && diffuseness < 1.00001f, "Erroneous parameter analysis");
335 doa_idx = pcon->doa_idx[band];
336 gain_idx = pcon->gains_idx[band];
337 gain_dir = pcon->gains_dir[band];
338 gain_diff = pcon->gains_diff[band];
339
340 /* Optional biasing (e.g. to conduct de-reverberation or to emphasise reverberation) */
341 streamBalance = SAF_CLAMP(s->streamBalance[band], 0.0f, 2.0f);
342 eq = s->eq[band];
343 if(streamBalance<1.0f){
344 a = streamBalance; /* pump more direct energy into output */
345 b = 1.0f; /* pass ambient stream as normal */
346 }
347 else {
348 a = 1.0f; /* pass source stream as normal */
349 b = 2.0f - streamBalance; /* pump less ambient energy into output */
350 }
351 a *= gain_dir;
352 b *= gain_diff;
353
354 /* Source array steering vector for the estimated DoAs */
355 for(i=0; i<nMics; i++)
356 s->As[i] = s->H_array[band*nMics*(s->nGrid) + i*(s->nGrid) + doa_idx];
357
358 /* Anechoic relative transfer functions (RTFs) */
359 for(i=0; i<nMics; i++){
360 s->As_l[i] = ccdivf(s->As[i], s->As[s->refIndices[0]]);
361 s->As_r[i] = ccdivf(s->As[i], s->As[s->refIndices[1]]);
362 }
363
364 /* HRTF for this reproduction DoA */
365 h_dir[0] = s->H_bin[band*NUM_EARS*(s->nGrid) + 0*(s->nGrid) + gain_idx];
366 h_dir[1] = s->H_bin[band*NUM_EARS*(s->nGrid) + 1*(s->nGrid) + gain_idx];
367 g_l = ccdivf(h_dir[0], s->As[s->refIndices[0]]); /* (Relative transfer functions) */
368 g_r = ccdivf(h_dir[1], s->As[s->refIndices[1]]);
369 if(cabsf(g_l)>4.0f || cabsf(g_r)>4.0f) /* if >12dB, then bypass: */
370 g_l = g_r = cmplxf(1.0f, 0.0f);
371
372 /* Diffuse mixing matrix (if the sound-field is analysed to be more diffuse, then we mix in more of just the reference sensors) */
373 memset(s->Q_diff, 0, NUM_EARS*nMics*sizeof(float_complex));
374 s->Q_diff[0*nMics+s->refIndices[0]] = cmplxf(s->diffEQ[band], 0.0f);
375 s->Q_diff[1*nMics+s->refIndices[1]] = cmplxf(s->diffEQ[band], 0.0f);
376
377 /* Source mixing matrix (beamforming towards the estimated DoAs) */
378 switch(s->beamOption){
379 case HADES_BEAMFORMER_NONE: /* No beamforming required */ break;
381 /* Normalise the beamformers to unity gain in the look direction */
382 utility_cpinv(s->hPinv, s->As_l, nMics, 1, s->Q_dir);
383 utility_cpinv(s->hPinv, s->As_r, nMics, 1, s->Q_dir + nMics);
384
385 /* Now bring their response from being w.r.t the array to being w.r.t the HRTF instead */
386 cblas_cscal(nMics, &g_l, s->Q_dir, 1);
387 cblas_cscal(nMics, &g_r, s->Q_dir + nMics, 1);
388 break;
389
391 /* prep */
392 cblas_ccopy(nMics*nMics, scon->Cx[band].Cx, 1, Cx, 1);
393 trace_M = 0.0f;
394 for(i=0; i<nMics; i++)
395 trace_M += crealf(Cx[i*nMics+i]);
396 sum_As = cblas_scasum(nMics, s->As, 1);
397
398 /* Compute beamforming weights if checks pass */
399 if( trace_M < 0.0001f || sum_As < 0.0001f)
400 memset(s->Q_dir, 0, NUM_EARS*nMics*sizeof(float_complex));
401 else{
402 /* Regularise Cx */
403 reg_M = (trace_M/(float)nMics) * 10.0f + 0.0001f;
404 for(i=0; i<nMics; i++)
405 Cx[i*nMics+i] = craddf(Cx[i*nMics+i], reg_M);
406
407 /* Compute MVDR weights w.r.t the reference sensor at each ear, [As^H Cx^-1 As]^-1 As^H Cx^-1 */
408 for(j=0; j<NUM_EARS; j++){
409 /* Solve As^H Cx-1 */
410 utility_cvconj(j==0 ? s->As_l : s->As_r, nMics, conj_As);
411 utility_cglslv(s->hLinSolve, Cx, nMics, conj_As, 1, AsH_invCx);
412
413 /* Compute As^H Cx-1 As */
414 utility_cvvdot(AsH_invCx, j==0 ? s->As_l : s->As_r, nMics, NO_CONJ, &AsH_invCx_As);
415 AsH_invCx_As = craddf(AsH_invCx_As, 0.00001f);
416
417 /* The solution */
418 AsH_invCx_As = ccdivf(cmplxf(1.0f, 0.0f), AsH_invCx_As);
419 cblas_cscal(nMics, &AsH_invCx_As, AsH_invCx, 1);
420 cblas_ccopy(nMics, AsH_invCx, 1, s->Q_dir + j*nMics, 1);
421 }
422
423 /* Now bring their response from being w.r.t the array to instead being w.r.t the HRTF */
424 cblas_cscal(nMics, &g_l, s->Q_dir, 1);
425 cblas_cscal(nMics, &g_r, s->Q_dir + nMics, 1);
426 }
427 break;
428 }
429
430 /* Prototype mixing matrix */
432 /* No beamforming (just pass through the reference signals) */
433 memset(s->Q, 0, NUM_EARS*nMics*sizeof(float_complex));
434 s->Q[0*nMics+s->refIndices[0]] = cmplxf(1.0f, 0.0f);
435 s->Q[1*nMics+s->refIndices[1]] = cmplxf(1.0f, 0.0f);
436 //s->Q[0*nMics+s->refIndices[0]] = cmplxf(s->diffEQ[band], 0.0f);
437 //s->Q[1*nMics+s->refIndices[1]] = cmplxf(s->diffEQ[band], 0.0f);
438 }
439 else{
440 /* Mix in the beamforming weights, conforming to the assumed direct-diffuse model */
441 cblas_ccopy(NUM_EARS*nMics, s->Q_dir, 1, s->Q, 1);
442 cblas_sscal(/*re+im*/2*NUM_EARS*nMics, eq*a*(1.0f-diffuseness), (float*)s->Q, 1);
443 cblas_saxpy(/*re+im*/2*NUM_EARS*nMics, eq*b*diffuseness, (float*)s->Q_diff, 1, (float*)s->Q, 1);
444 }
445
446 /* Target output signal energy (used for the covariance matching) */
447 targetEnergy = 0.0f;
448 for(i=0; i<nMics; i++)
449 targetEnergy += crealf(scon->Cx[band].Cx[i*nMics+i]);
450 targetEnergy = eq*0.25f*targetEnergy * s->diffEQ[band];
451
452 /* Final mixing matrix */
453 if(s->enableCM && targetEnergy>0.0001f){
454 /* "Direct" contributions to the target spatial covariance matrix */
455 cblas_cgemm(CblasRowMajor, CblasNoTrans, CblasConjTrans, NUM_EARS, NUM_EARS, 1, &calpha,
456 h_dir, 1,
457 h_dir, 1, &cbeta,
458 s->Cy, NUM_EARS);
459 cblas_sscal(/*re+im*/2*NUM_EARS*NUM_EARS, eq*a*(1.0f-diffuseness)*targetEnergy, (float*)s->Cy, 1);
460
461 /* "Diffuse" contributions to the target spatial covariance matrix */
462 cblas_saxpy(/*re+im*/2*NUM_EARS*NUM_EARS, eq*b*diffuseness*targetEnergy, (float*)&(s->DCM_bin_norm[band*NUM_EARS*NUM_EARS]), 1, (float*)s->Cy, 1);
463
464 /* Solve the covariance matching problem */
465 formulate_M_and_Cr_cmplx(s->hCDF, (float_complex*)scon->Cx[band].Cx, s->Cy, s->Q, 1, 0.1f, s->new_M, NULL);
466 }
467 else
468 cblas_ccopy(NUM_EARS*nMics, s->Q, 1, s->new_M, 1);
469
470 /* Optional Equalisation */
471 cblas_sscal(/*re+im*/2*NUM_EARS*nMics, eq, (float*)s->new_M, 1);
472
473 /* Temporal averaging of mixing matrices */
474 cblas_sscal(/*re+im*/2*NUM_EARS*nMics, synAvgCoeff, (float*)s->M[band], 1);
475 cblas_saxpy(/*re+im*/2*NUM_EARS*nMics, 1.0f-synAvgCoeff, (float*)s->new_M, 1, (float*)s->M[band], 1);
476 }
477
478 /* Apply mixing matrices */
479 for(band=0; band<s->nBands; band++){
480 cblas_cgemm(CblasRowMajor, CblasNoTrans, CblasNoTrans, NUM_EARS, s->timeSlots, nMics, &calpha,
481 s->M[band], nMics,
482 FLATTEN2D(scon->inTF[band]), s->timeSlots, &cbeta,
483 FLATTEN2D(s->outTF[band]), s->timeSlots);
484 }
485
486 /* inverse time-frequency transform */
487 switch(s->fbOpt){
488 case HADES_USE_AFSTFT_LD: /* fall through */
490 }
491
492 /* Copy to output */
493 for(ch=0; ch<SAF_MIN(nChannels, NUM_EARS); ch++)
494 memcpy(output[ch], s->outTD[ch], blocksize*sizeof(float));
495 for(; ch<nChannels; ch++)
496 memset(output[ch], 0, blocksize*sizeof(float));
497}
498
500(
501 hades_synthesis_handle const hSyn,
502 int* nBands
503)
504{
506 if(hSyn==NULL){
507 if(nBands!=NULL)
508 (*nBands) = 0;
509 return NULL;
510 }
511 s = (hades_synthesis_data*)(hSyn);
512 if(nBands!=NULL)
513 (*nBands) = s->nBands;
514 return s->eq;
515}
516
518(
519 hades_synthesis_handle const hSyn,
520 int* nBands
521)
522{
524 if(hSyn==NULL){
525 if(nBands!=NULL)
526 (*nBands) = 0;
527 return NULL;
528 }
529 s = (hades_synthesis_data*)(hSyn);
530 if(nBands!=NULL)
531 (*nBands) = s->nBands;
532 return s->streamBalance;
533}
534
536(
537 hades_synthesis_handle const hSyn
538)
539{
541 if(hSyn==NULL)
542 return NULL;
543 s = (hades_synthesis_data*)(hSyn);
544 return &(s->synAvgCoeff);
545}
547(
548 hades_synthesis_handle const hSyn
549)
550{
551 if(hSyn==NULL)
552 return 0;
553 return 0; /* Accounted for in hades_analysis_getProcDelay() */
554}
555
556#endif /* SAF_ENABLE_HADES_MODULE */
void afSTFT_clearBuffers(void *const hSTFT)
Flushes time-domain buffers with zeros.
Definition afSTFTlib.c:519
void afSTFT_create(void **const phSTFT, int nCHin, int nCHout, int hopsize, int lowDelayMode, int hybridmode, AFSTFT_FDDATA_FORMAT format)
Creates an instance of afSTFT.
Definition afSTFTlib.c:143
void afSTFT_backward_knownDimensions(void *const hSTFT, float_complex ***dataFD, int framesize, int dataFD_nCH, int dataFD_nHops, float **dataTD)
Performs backward afSTFT transform (dataFD dimensions are known)
Definition afSTFTlib.c:391
void afSTFT_destroy(void **const phSTFT)
Destroys an instance of afSTFT.
Definition afSTFTlib.c:199
@ AFSTFT_BANDS_CH_TIME
nBands x nChannels x nTimeHops
Definition afSTFTlib.h:80
void formulate_M_and_Cr_cmplx(void *const hCdf, float_complex *Cx, float_complex *Cy, float_complex *Q, int useEnergyFLAG, float reg, float_complex *M, float_complex *Cr)
Computes the optimal mixing matrices.
void cdf4sap_cmplx_destroy(void **const phCdf)
Destroys an instance of the Covariance Domain Framework.
void cdf4sap_cmplx_create(void **const phCdf, int nXcols, int nYcols)
Creates an instance of the Covariance Domain Framework.
#define saf_assert(x, message)
Macro to make an assertion, along with a string explaining its purpose.
void utility_cpinv_create(void **const phWork, int maxDim1, int maxDim2)
(Optional) Pre-allocate the working struct used by utility_cpinv()
#define SAF_CLAMP(a, min, max)
Ensures value "a" is clamped between the "min" and "max" values.
void utility_cglslv_destroy(void **const phWork)
De-allocate the working struct used by utility_cglslv()
void utility_cglslv(void *const hWork, const float_complex *A, const int dim, float_complex *B, int nCol, float_complex *X)
General linear solver: single precision complex, i.e.
#define NUM_EARS
2 (true for most humans)
void utility_cpinv(void *const hWork, const float_complex *inM, const int dim1, const int dim2, float_complex *outM)
General matrix pseudo-inverse (the svd way): single precision complex, i.e.
void utility_cpinv_destroy(void **const phWork)
De-allocate the working struct used by utility_cpinv()
#define SAF_MIN(a, b)
Returns the minimum of the two values.
void utility_cvvdot(const float_complex *a, const float_complex *b, const int len, CONJ_FLAG flag, float_complex *c)
Single-precision, complex, vector-vector dot product, i.e.
void utility_cglslv_create(void **const phWork, int maxDim, int maxNCol)
(Optional) Pre-allocate the working struct used by utility_cglslv()
void utility_cvconj(const float_complex *a, const int len, float_complex *c)
Single-precision, complex, vector-conjugate, i.e.
@ NO_CONJ
Do not take the conjugate.
void ** malloc2d(size_t dim1, size_t dim2, size_t data_size)
2-D malloc (contiguously allocated, so use free() as usual to deallocate)
Definition md_malloc.c:89
void * malloc1d(size_t dim1_data_size)
1-D malloc (same as malloc, but with error checking)
Definition md_malloc.c:59
void * calloc1d(size_t dim1, size_t data_size)
1-D calloc (same as calloc, but with error checking)
Definition md_malloc.c:69
void *** malloc3d(size_t dim1, size_t dim2, size_t dim3, size_t data_size)
3-D malloc (contiguously allocated, so use free() as usual to deallocate)
Definition md_malloc.c:151
#define FLATTEN2D(A)
Use this macro when passing a 2-D dynamic multi-dimensional array to memset, memcpy or any other func...
Definition md_malloc.h:65
#define HADES_MAX_NMICS
Maximum number of microphones.
struct _hades_analysis_data * hades_analysis_handle
Handle for the hades analysis data.
struct _hades_param_container_data * hades_param_container_handle
Handle for the hades parameter container data.
struct _hades_signal_container_data * hades_signal_container_handle
Handle for the hades signal container data.
@ HADES_USE_AFSTFT
Alias-free STFT filterbank.
@ HADES_USE_AFSTFT_LD
Alias-free STFT filterbank (low delay)
void hades_getInterpolatedHRTFs(hades_analysis_handle const hAna, HADES_HRTF_INTERP_OPTIONS interpOption, hades_binaural_config *binConfig, float *target_dirs_deg, int nTargetDirs, float_complex *hrtf_interp)
Binaural filter interpolator.
Internal header for the HADES module (SAF_HADES_MODULE)
void hades_radial_editor_destroy(hades_radial_editor_handle *const phREd)
Destroys an instance of a hades radial editor object.
int hades_synthesis_getProcDelay(hades_synthesis_handle const hSyn)
Returns the synthesiser processing delay, in samples.
void hades_radial_editor_apply(hades_radial_editor_handle const hREd, hades_param_container_handle const hPCon, float dirGain_dB[360])
Applies the radial (360 degree) parameter editing.
void hades_radial_editor_create(hades_radial_editor_handle *const phREd, hades_analysis_handle const hAna)
Creates and returns a handle to an instance of a hades radial editor object, which allows for directi...
float * hades_synthesis_getEqPtr(hades_synthesis_handle const hSyn, int *nBands)
Returns a pointer to the eq vector, which can be changed at run-time.
float * hades_synthesis_getSynthesisAveragingCoeffPtr(hades_synthesis_handle const hSyn)
Returns a pointer to the synthesis averaging coefficient scalar [0..1], which can be changed at run-t...
void hades_synthesis_apply(hades_synthesis_handle const hSyn, hades_param_container_handle const hPCon, hades_signal_container_handle const hSCon, int nChannels, int blocksize, float **output)
Performs hades synthesis.
void hades_synthesis_create(hades_synthesis_handle *const phSyn, hades_analysis_handle const hAna, HADES_BEAMFORMER_TYPE beamOption, int enableCM, int refIndices[2], hades_binaural_config *binConfig, HADES_HRTF_INTERP_OPTIONS interpOption)
Creates and returns a handle to an instance of a hades synthesis object.
float * hades_synthesis_getStreamBalancePtr(hades_synthesis_handle const hSyn, int *nBands)
Returns a pointer to the stream balance vector [0..2], which can be changed at run-time.
void hades_synthesis_reset(hades_synthesis_handle const hSyn)
Flushes run-time buffers with zeros.
void hades_synthesis_destroy(hades_synthesis_handle *const phSyn)
Destroys an instance of hades synthesis.
Header for the HADES synthesis (SAF_HADES_MODULE)
struct _hades_radial_editor_data * hades_radial_editor_handle
Handle for the hades radial editor data.
struct _hades_synthesis_data * hades_synthesis_handle
Handle for the hades synthesis data.
HADES_BEAMFORMER_TYPE
Beamforming options for hades_synthesis.
@ HADES_BEAMFORMER_FILTER_AND_SUM
Filter-and-sum beamforming.
@ HADES_BEAMFORMER_BMVDR
Binaural minimum-variance distortion- less response (MVDR) beamforming.
@ HADES_BEAMFORMER_NONE
No beamforming (ref sensors only)
HADES_HRTF_INTERP_OPTIONS
HRTF interpolation options for hades_synthesis.
Main structure for hades analysis.
int nMics
Number of microphones.
int blocksize
Number of samples to process at a time (note that 1 doa and diffuseness estimate is made per block)
int hybridmode
Optionally, the lowest TF bands may be subdivided to improve low-freq resolution.
float_complex * DCM_array
Diffuse covariance matrix (computed over all grid directions and weighted); FLAT: nBands x nMics x nM...
float * grid_dirs_xyz
Scanning grid coordinates (unit vectors and only used by grid-based estimators); FLAT: nGrid x 3.
float * freqVector
Centre frequencies; nBands x 1.
float * grid_dirs_deg
Array grid dirs in degrees; FLAT: nGrid x 2.
int hopsize
Filterbank hop size (blocksize must be divisable by this.
HADES_FILTERBANKS fbOpt
see HADES_FILTERBANKS
float_complex * W
Diffuse integration weighting matrix; FLAT: nGrid x nGrid.
float_complex * H_array
Array IRs in the frequency domain; FLAT: nBands x nMics x nDirs.
int nBands
Number of frequency bands.
int timeSlots
Number of time slots.
int nGrid
Number of grid/scanning directions.
Binaural configuration struct.
int nHRIR
Number of HRIRs.
int hrir_fs
HRIR sample rate.
int lHRIR
Length of HRIRs in samples.
float * hrir_dirs_deg
HRTF directions in [azimuth elevation] format, in degrees; FLAT: nHRIR x 2.
float * hrirs
Matrix of HRIR data; FLAT: nHRIR x NUM_EARS x lHRIR.
Parameter container to store the data from an analyser for one blocksize of audio.
int * gains_idx
Reproduction direction index per band; nBands x 1.
int * doa_idx
Beamforming direction index per band; nBands x 1.
float * diffuseness
Diffuseness value per band; nBands x 1.
float * gains_diff
Extra diffuse reproduction gain per band (default=1.0f); nBands x 1
float * gains_dir
Extra direct reproduction gain per band (default=1.0f); nBands x 1
Main structure for hades radial (360degree) gain and direct-to-diffuse ratio editor.
float * pGrid_dirs_deg
Pointer to grid dirs in degrees; FLAT: nGrid x 2.
float * pGrid_dirs_xyz
Pointer to grid dirs as Cartesian coordinates of unit length; FLAT: nGrid x 3.
int nGrid
Number of grid/scanning directions.
Signal container to store one block of TF-domain audio data.
CxMic * Cx
NON-time-averaged covariance matrix per band; nBands x .Cx(nMics x nMics)
float_complex *** inTF
Input frame in TF-domain; nBands x nMics x timeSlots.
Main structure for hades synthesis.
float ** outTD
output time-domain buffer; NUM_EARS x blocksize
void * hLinSolve
Handle for solving linear equations (Ax=b)
int blocksize
blocksize in samples
void * hCDF
Handle for solving the covariance matching problem.
float_complex ** M
Mixing matrix per band; nBands x FLAT: (NUM_EARS x nMics)
float synAvgCoeff
Mixing matrix averaging coefficent [0..1].
HADES_FILTERBANKS fbOpt
Filterbank option, see HADES_FILTERBANKS.
float_complex * new_M
New mixing matrix (not yet temporally averaged); FLAT: NUM_EARS x nMics.
int nBands
Number of bands in the time-frequency transform domain.
float_complex * H_array
Array IRs in the frequency domain; FLAT: nBands x nMics x nGrid.
float_complex * As
Array steering vector for DoA; FLAT: nMics x 1.
float * eq
Gain factor per band; nBands x 1.
int timeSlots
Number of time frames in the time-frequency transform domain.
float_complex * H_bin
To spatialise the source beamformers; FLAT: nBands x NUM_EARS x nGrid.
HADES_BEAMFORMER_TYPE beamOption
see HADES_BEAMFORMER_TYPE
int enableCM
Flag: whether the spatial covariance matching is enabled (1) or disabled (0)
hades_binaural_config * binConfig
Internal copy of user configuration.
float * freqVector
Frequency vector (band centre frequencies); nBands x 1.
float * grid_dirs_deg
Array grid dirs in degrees; FLAT: nGrid x 2.
int nGrid
Number of grid/scanning directions.
float_complex * Q
Mixing matrix for the direct and diffuse streams combined (based on the diffuseness value); FLAT: NUM...
float_complex * As_r
Array steering vector relative to right reference sensor; FLAT: nMics x 1.
float * diffEQ
EQ curve to bring the overall diffuse-field magnitude response of the array to that of the HRTFs inst...
void * hPinv
Handle for computing the Moore-Penrose pseudo inverse.
float_complex * DCM_array
Diffuse coherence matrix for the array; FLAT: nBands x nMics x nMics.
int hopsize
hopsize in samples
float_complex * Cy
Target binaural spatial covariance matrix; FLAT: NUM_EARS x NUM_EARS.
float_complex *** outTF
nBands x NUM_EARS x timeSlots
float ** grid_dirs_xyz
Grid dirs as Cartesian coordinates of unit length; nGrid x 3.
float_complex * Q_dir
Mixing matrix for the direct stream; FLAT: NUM_EARS x nMics.
float_complex * DCM_bin_norm
Diffuse coherence matrix for the HRTF set, normalised with 1/trace(DCM_bin); FLAT: nBands x nMics x n...
float_complex * W
Diffuse integration weighting matrix; FLAT: nGrid x nGrid.
float_complex * Q_diff
Mixing matrix for the diffuse stream; FLAT: NUM_EARS x nMics.
HADES_HRTF_INTERP_OPTIONS interpOption
HRIR interpolation option, see HADES_HRTF_INTERP_OPTIONS.
float_complex * As_l
Array steering vector relative to left reference sensor; FLAT: nMics x 1.
int refIndices[2]
Indices into [0 nMics-1], defining the reference sensors.
void * hFB_dec
Filterbank handle.
int nMics
Number of microphones.
float * streamBalance
Stream balance per band (0:fully diffuse, 1:balanced, 2:fully direct); nBands x 1.