[6db11f]: mex / comp_ifilterbank_fftbl.c  Maximize  Restore  History

Download this file

176 lines (134 with data), 4.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
#ifndef _LTFAT_MEX_FILE
#define _LTFAT_MEX_FILE
#define ISNARGINEQ 5
#define TYPEDEPARGS 0, 1
#define SINGLEARGS
#define COMPLEXARGS
#define EXPORTALIAS comp_ifilterbank_fftbl
#endif // _LTFAT_MEX_FILE - INCLUDED ONCE
#define MEX_FILE __BASE_FILE__
#include "ltfat_mex_template_helper.h"
#if defined(LTFAT_SINGLE) || defined(LTFAT_DOUBLE)
#include "ltfat_types.h"
#include <math.h>
#include "config.h"
static LTFAT_FFTW(plan)** LTFAT_NAME(oldPlans) = 0;
static mwSize* LTFAT_NAME(oldLc) = 0;
static mwSize LTFAT_NAME(oldM) = 0;
void LTFAT_NAME(fftblMexAtExitFnc)()
{
if(LTFAT_NAME(oldPlans)!=0)
{
for(mwIndex m=0;m<LTFAT_NAME(oldM);m++)
{
if(LTFAT_NAME(oldPlans)[m]!=0)
{
LTFAT_FFTW(destroy_plan)(*LTFAT_NAME(oldPlans)[m]);
ltfat_free(LTFAT_NAME(oldPlans)[m]);
}
}
ltfat_free(LTFAT_NAME(oldPlans));
LTFAT_NAME(oldPlans) = 0;
}
if(LTFAT_NAME(oldLc)!=0)
{
ltfat_free(LTFAT_NAME(oldLc));
LTFAT_NAME(oldLc) = 0;
}
LTFAT_NAME(oldM) = 0;
#ifdef _DEBUG
mexPrintf("Exit fnc called: %s\n",__PRETTY_FUNCTION__);
#endif
}
// Calling convention:
// c = comp_ifilterbank_fftbl(c,G,foff,a,realonly)
void LTFAT_NAME(ltfatMexFnc)( int nlhs, mxArray *plhs[],int nrhs, const mxArray *prhs[] )
{
static int atExitFncRegistered = 0;
if(!atExitFncRegistered)
{
LTFAT_NAME(ltfatMexAtExit)(LTFAT_NAME(fftblMexAtExitFnc));
atExitFncRegistered = 1;
}
const mxArray* mxc = prhs[0];
const mxArray* mxG = prhs[1];
double* foff = (double*) mxGetData(prhs[2]);
double* a = (double*) mxGetData(prhs[3]);
double* realonly = (double*) mxGetData(prhs[4]);
// number of channels
mwSize W = mxGetN(mxGetCell(mxc,0));
// filter number
mwSize M = mxGetNumberOfElements(mxG);
//
mwSize acols = mxGetN(prhs[3]);
double afrac[M];
memcpy(afrac,a,M*sizeof(double));
if(acols>1)
{
for(mwIndex m=0;m<M;m++)
{
afrac[m] = afrac[m]/a[M+m];
}
}
// POINTER TO THE FILTERS
LTFAT_REAL _Complex* GPtrs[M];
// input lengths
mwSize inLen[M];
// POINTER TO INPUTS
LTFAT_REAL _Complex* cPtrs[M];
// filter lengths
mwSize filtLen[M];
for(mwIndex m = 0; m < M; m++)
{
cPtrs[m] = (LTFAT_REAL _Complex*) mxGetData(mxGetCell(mxc,m));
inLen[m] = (mwSize) mxGetM(mxGetCell(mxc, m));
GPtrs[m] = (LTFAT_REAL _Complex*) mxGetData(mxGetCell(mxG, m));
LTFAT_NAME_COMPLEX(conjugate_array);
filtLen[m] = (mwSize) mxGetNumberOfElements(mxGetCell(mxG, m));
}
// output data length
mwSize L = (mwSize) floor(afrac[0]*inLen[0] + 0.5);
/*
WOW. You really check for that?
*/
// if(plhs[0]==NULL || ~mxIsNumeric(plhs[0]))
// {
plhs[0] = ltfatCreateMatrix(L, W, LTFAT_MX_CLASSID, mxCOMPLEX);
memset(mxGetData(plhs[0]),0,L*W*sizeof(LTFAT_REAL _Complex*));
// }
mxArray* mxF = plhs[0];
LTFAT_REAL _Complex* FPtr = (LTFAT_REAL _Complex*) mxGetData(mxF);
if(M!=LTFAT_NAME(oldM))
{
LTFAT_NAME(fftblMexAtExitFnc)();
LTFAT_NAME(oldM) = M;
LTFAT_NAME(oldLc) = (mwSize*) ltfat_calloc(M,sizeof(mwSize));
LTFAT_NAME(oldPlans) = (LTFAT_FFTW(plan)**) ltfat_calloc(M,sizeof(LTFAT_FFTW(plan)*));
}
// over all channels
for(mwIndex m =0; m<M; m++)
{
LTFAT_REAL _Complex* buffer = (LTFAT_REAL _Complex*) ltfat_malloc(inLen[m]*sizeof(LTFAT_REAL _Complex));
if(LTFAT_NAME(oldLc)[m]!=inLen[m])
{
LTFAT_NAME(oldLc)[m] = inLen[m];
LTFAT_FFTW(plan) ptmp = LTFAT_FFTW(plan_dft_1d)(inLen[m],(LTFAT_REAL (*)[2]) buffer,(LTFAT_REAL (*)[2]) buffer, FFTW_FORWARD, FFTW_OPTITYPE);
if(LTFAT_NAME(oldPlans)[m]!=0)
{
LTFAT_FFTW(destroy_plan)(*LTFAT_NAME(oldPlans)[m]);
ltfat_free(LTFAT_NAME(oldPlans)[m]);
}
LTFAT_NAME(oldPlans)[m] = ltfat_malloc(sizeof(ptmp));
memcpy(LTFAT_NAME(oldPlans)[m],&ptmp,sizeof(ptmp));
}
for(mwIndex w =0; w<W; w++)
{
// Obtain pointer to w-th column in output and input
LTFAT_REAL _Complex *FPtrCol = FPtr + w*L;
LTFAT_REAL _Complex *cPtrCol = cPtrs[m] + w*inLen[m];
LTFAT_NAME(upconv_fftbl_plan)(cPtrCol,inLen[m],GPtrs[m],filtLen[m],(int)foff[m],afrac[m],realonly[m],FPtrCol,LTFAT_NAME(oldPlans)[m],buffer);
}
ltfat_free(buffer);
}
}
#endif

Get latest updates about Open Source Projects, Conferences and News.

Sign up for the SourceForge newsletter:





No, thanks