Learn how easy it is to sync an existing GitHub or Google Code repo to a SourceForge project! See Demo

Close

[e478b6]: mex / comp_filterbank.c Maximize Restore History

Download this file

comp_filterbank.c    366 lines (299 with data), 10.3 kB

  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
#include "mex.h"
#include "fftw3.h"
#include <string.h>
void comp_filterbank_td(int nlhs, mxArray *plhs[],int nrhs, const mxArray *prhs[] );
void comp_filterbank_fft(int nlhs, mxArray *plhs[],int nrhs, const mxArray *prhs[] );
void comp_filterbank_fftbl(int nlhs, mxArray *plhs[],int nrhs, const mxArray *prhs[] );
void comp_filterbank_fftbl_atexit();
void comp_filterbank_fft_atexit();
void comp_filterbank_td_atexit();
static fftw_plan* p_double = NULL;
static fftwf_plan* p_float = NULL;
/*
Since the array is store for the lifetime of the MEX, we introduce limit od the array length.
2^20 ~ 16 MB of complex double
*/
#define MAXARRAYLEN 1048576
// Static pointer for holding the array the FFTW plan is
static mxArray* mxF = NULL;
// Calling convention:
// comp_filterbank(f,g,a);
/*
MEX exit fnc. are not called by Matlab
*/
void filterbankAtExit()
{
#ifdef _DEBUG
mexPrintf("Exit fnc called: %s\n",__PRETTY_FUNCTION__);
#endif
comp_filterbank_fftbl_atexit();
comp_filterbank_fft_atexit();
comp_filterbank_td_atexit();
if(mxF!=NULL)
mxDestroyArray(mxF);
if(p_double!=NULL)
{
fftw_destroy_plan(*p_double);
free(p_double);
}
if(p_float!=NULL)
{
fftw_destroy_plan(*p_float);
free(p_float);
}
}
void mexFunction( int nlhs, mxArray *plhs[],int nrhs, const mxArray *prhs[] )
{
const mxArray* mxf = prhs[0];
const mxArray* mxg = prhs[1];
const mxArray* mxa = prhs[2];
// input data length
const mwSize L = mxGetM(mxf);
const mwSize W = mxGetN(mxf);
// filter number
const mwSize M = mxGetNumberOfElements(mxg);
// a col count
mwSize acols = mxGetN(mxa);
// pointer to a
double *a = (double*) mxGetData(mxa);
if(acols>1)
{
int isOnes = 1;
for(mwIndex m=0;m<M;m++)
{
isOnes = isOnes && a[M+m]==1;
}
if(isOnes)
{
acols = 1;
}
}
// Cell output
plhs[0] = mxCreateCellMatrix(M, 1);
// Stuff for sorting the filters
mwSize tdCount = 0;
mwSize fftCount = 0;
mwSize fftblCount = 0;
mwIndex tdArgsIdx[M];
mwIndex fftArgsIdx[M];
mwIndex fftblArgsIdx[M];
// WALK the filters to determine what has to be done
for(mwIndex m =0; m<M; m++)
{
mxArray * gEl = mxGetCell(mxg, m);
if(mxGetField(gEl,0,"h")!=NULL)
{
tdArgsIdx[tdCount++] = m;
continue;
}
if(mxGetField(gEl,0,"H")!=NULL)
{
if(acols==1&&L==mxGetNumberOfElements(mxGetField(gEl,0,"H")))
{
fftArgsIdx[fftCount++] = m;
continue;
}
else
{
fftblArgsIdx[fftblCount++] = m;
continue;
}
}
}
if(tdCount>0)
{
/*
Here, we have to reformat the inputs and pick up results to comply with:
c=comp_filterbank_td(f,g,a,offset,ext);
BEWARE OF THE AUTOMATIC DEALLOCATION!! by the Matlab engine.
Arrays can be very easily freed twice causing segfaults.
This happends particulary when using mxCreateCell* which stores
pointers to other mxArray structs. Setting all such pointers to
NULL after they are used seems to solve it.
*/
mxArray* plhs_td[1];
const mxArray* prhs_td[5];
prhs_td[0] = mxf;
prhs_td[1] = mxCreateCellMatrix(tdCount,1);
prhs_td[2] = mxCreateDoubleMatrix(tdCount,1,mxREAL);
double* aPtr = (double*)mxGetPr(prhs_td[2]);
prhs_td[3] = mxCreateDoubleMatrix(tdCount,1,mxREAL);
double* offsetPtr = (double*)mxGetPr(prhs_td[3]);
prhs_td[4] = mxCreateString("per");
for(mwIndex m=0;m<tdCount;m++)
{
mxArray * gEl = mxGetCell(mxg, tdArgsIdx[m]);
mxSetCell((mxArray*)prhs_td[1],m,mxGetField(gEl,0,"h"));
// This has overhead
//mxSetCell((mxArray*)prhs_td[1],m,mxDuplicateArray(mxGetField(gEl,0,"h")));
aPtr[m] = a[tdArgsIdx[m]];
offsetPtr[m] = mxGetScalar(mxGetField(gEl,0,"offset"));
}
// Finally call it!
comp_filterbank_td(1,plhs_td,5, prhs_td);
// This has overhead:
// mexCallMATLAB(1,plhs_td,5, prhs_td,"comp_filterbank_td");
// Copy pointers to a proper index in the output + unset all duplicate cell elements
for(mwIndex m=0;m<tdCount;m++)
{
mxSetCell(plhs[0],tdArgsIdx[m],mxGetCell(plhs_td[0],m));
mxSetCell(plhs_td[0],m,NULL);
mxSetCell((mxArray*)prhs_td[1],m,NULL);
}
mxDestroyArray((mxArray*)plhs_td[0]);
mxDestroyArray((mxArray*)prhs_td[1]);
mxDestroyArray((mxArray*)prhs_td[2]);
mxDestroyArray((mxArray*)prhs_td[3]);
mxDestroyArray((mxArray*)prhs_td[4]);
}
if(fftCount>0 || fftblCount>0)
{
// Need to do FFT of mxf
mwIndex ndim = 2;
const mwSize dims[] = {L,W};
if(mxF==NULL || mxGetM(mxF)!=L || mxGetN(mxF)!=W)
{
if(mxF!=NULL)
{
mxDestroyArray(mxF);
mxF = NULL;
printf("Should be called just once\n");
}
if(mxIsDouble(mxf))
{
mxF = mxCreateNumericArray(ndim,dims,mxDOUBLE_CLASS,mxCOMPLEX);
fftw_iodim fftw_dims[1];
fftw_iodim howmanydims[1];
fftw_dims[0].n = L;
fftw_dims[0].is = 1;
fftw_dims[0].os = 1;
howmanydims[0].n = W;
howmanydims[0].is = L;
howmanydims[0].os = L;
if(p_double==NULL)
p_double = (fftw_plan*) malloc(sizeof(fftw_plan));
else
fftw_destroy_plan(*p_double);
*p_double = fftw_plan_guru_split_dft(
1, fftw_dims,
1, howmanydims,
mxGetPr(mxF), mxGetPi(mxF), mxGetPr(mxF), mxGetPi(mxF),
FFTW_MEASURE);
}
else if(mxIsSingle(mxf))
{
mxF = mxCreateNumericArray(ndim,dims,mxSINGLE_CLASS,mxCOMPLEX);
// mexPrintf("M= %i, N= %i\n",mxGetM(mxF),mxGetN(mxF));
fftwf_iodim fftw_dims[1];
fftwf_iodim howmanydims[1];
fftw_dims[0].n = L;
fftw_dims[0].is = 1;
fftw_dims[0].os = 1;
howmanydims[0].n = W;
howmanydims[0].is = L;
howmanydims[0].os = L;
if(p_float==NULL)
p_float = (fftwf_plan*) malloc(sizeof(fftwf_plan));
else
fftwf_destroy_plan(*p_float);
*p_float = fftwf_plan_guru_split_dft(
1, fftw_dims,
1, howmanydims,
(float*)mxGetPr(mxF), (float*)mxGetPi(mxF),
(float*) mxGetPr(mxF), (float*)mxGetPi(mxF),
FFTW_ESTIMATE);
}
}
if(mxIsDouble(mxf))
{
memcpy(mxGetPr(mxF),mxGetPr(mxf),L*W*sizeof(double));
memset(mxGetPi(mxF),0,L*W*sizeof(double));
if(mxIsComplex(mxf))
memcpy(mxGetPi(mxF),mxGetPi(mxf),L*W*sizeof(double));
fftw_execute(*p_double);
}
else if(mxIsSingle(mxf))
{
memcpy(mxGetPr(mxF),mxGetPr(mxf),L*W*sizeof(float));
memset(mxGetPi(mxF),0,L*W*sizeof(float));
if(mxIsComplex(mxf))
memcpy(mxGetPi(mxF),mxGetPi(mxf),L*W*sizeof(float));
fftwf_execute(*p_float);
}
}
if(fftCount>0)
{
mxArray* plhs_fft[1];
const mxArray* prhs_fft[3];
prhs_fft[0] = mxF;
prhs_fft[1] = mxCreateCellMatrix(fftCount,1);
prhs_fft[2] = mxCreateDoubleMatrix(fftCount,1,mxREAL);
double* aPtr = (double*)mxGetPr(prhs_fft[2]);
for(mwIndex m=0;m<fftCount;m++)
{
mxArray * gEl = mxGetCell(mxg, fftArgsIdx[m]);
mxSetCell((mxArray*)prhs_fft[1],m,mxGetField(gEl,0,"H"));
// This has overhead
//mxSetCell((mxArray*)prhs_td[1],m,mxDuplicateArray(mxGetField(gEl,0,"h")));
aPtr[m] = a[fftArgsIdx[m]];
}
comp_filterbank_fft(1,plhs_fft,3, prhs_fft);
for(mwIndex m=0;m<fftCount;m++)
{
mxSetCell(plhs[0],fftArgsIdx[m],mxGetCell(plhs_fft[0],m));
mxSetCell(plhs_fft[0],m,NULL);
mxSetCell((mxArray*)prhs_fft[1],m,NULL);
}
mxDestroyArray((mxArray*)plhs_fft[0]);
mxDestroyArray((mxArray*)prhs_fft[1]);
mxDestroyArray((mxArray*)prhs_fft[2]);
}
if(fftblCount>0)
{
mxArray* plhs_fftbl[1];
const mxArray* prhs_fftbl[5];
prhs_fftbl[0] = mxF;
prhs_fftbl[1] = mxCreateCellMatrix(fftblCount,1);
prhs_fftbl[2] = mxCreateDoubleMatrix(fftblCount,1,mxREAL);
prhs_fftbl[3] = mxCreateDoubleMatrix(fftblCount,2,mxREAL);
prhs_fftbl[4] = mxCreateDoubleMatrix(fftblCount,1,mxREAL);
double* foffPtr = (double*)mxGetPr(prhs_fftbl[2]);
double* aPtr = (double*)mxGetPr(prhs_fftbl[3]);
double* realonlyPtr = (double*)mxGetPr(prhs_fftbl[4]);
for(mwIndex m=0;m<fftblCount;m++)
{
mxArray * gEl = mxGetCell(mxg, fftblArgsIdx[m]);
mxSetCell((mxArray*)prhs_fftbl[1],m,mxGetField(gEl,0,"H"));
foffPtr[m] = mxGetScalar(mxGetField(gEl,0,"foff"));
aPtr[m] = a[fftblArgsIdx[m]];
if(acols>1)
aPtr[m+fftblCount] = a[fftblArgsIdx[m]+M];
else
aPtr[m+fftblCount] = 1;
realonlyPtr[m] = mxGetScalar(mxGetField(gEl,0,"realonly"));
}
comp_filterbank_fftbl(1,plhs_fftbl,5, prhs_fftbl);
for(mwIndex m=0;m<fftblCount;m++)
{
mxSetCell(plhs[0],fftblArgsIdx[m],mxGetCell(plhs_fftbl[0],m));
mxSetCell(plhs_fftbl[0],m,NULL);
mxSetCell((mxArray*)prhs_fftbl[1],m,NULL);
}
mxDestroyArray((mxArray*)plhs_fftbl[0]);
mxDestroyArray((mxArray*)prhs_fftbl[1]);
mxDestroyArray((mxArray*)prhs_fftbl[2]);
mxDestroyArray((mxArray*)prhs_fftbl[3]);
mxDestroyArray((mxArray*)prhs_fftbl[4]);
}
/* This should overwrite function registered by mexAtExit in any of the previously
called MEX files */
mexAtExit(filterbankAtExit);
mexMakeArrayPersistent(mxF);
if(L*W>MAXARRAYLEN)
{
//printf("Damn. Should not get here\n");
mxDestroyArray(mxF);
mxF = NULL;
}
//int prd = 0;
}