SigPack - the C++ signal processing library
 
Loading...
Searching...
No Matches
Matlab

General

Armadillo provides a wrapper for Matlab integration in .../mex_interface/armaMex.hpp. Your function must first be compiled to a mex file by the Matlab compiler.

The way forward is :

  • create your c++ wrapper e.g myFunc.cpp (the file name will also be the function name in Matlab)
  • compile it using mex in Matlab
  • put the mex file where Matlab can find it

Compilation

Compilation is done in the Matlab environment. Compile parameters/flags is basically the same as when you build the "normal" SigPack applications, so if you want FFTW support you need to add e.g -DHAVE_FFTW -lfftw3

Windows

mex -I<Armadillo install dir>\include ...
-I<Armadillo install dir>\mex_interface ...
-I<Sigpack install dir>\include ...
-L<Armadillo install dir>\examples\lib_win64 ...
-lblas_win64_MT.lib -llapack_win64_MT.lib ...
myFunc.cpp

Linux

mex -I<Armadillo install dir>/include ...
-I<Armadillo install dir>/mex_interface ...
-I<Sigpack install dir>/include ...
-larmadillo ...
myFunc.cpp

Matlab code

% Demo
X = 0.5+randn(10000,1); % Generate random signal
Y = myFunc(X); % Run the demo using your c++ code
% Plot spectrum
figure;
plot(10*log10(abs(Y)),'b')
Y2 = pwelch(filter(fir1(17,0.25),1,X),512,256,'twosided','power');
hold on
plot(10*log10(abs(Y2)),'+r')
legend('SigPack','Matlab')

C++ code

#include "armaMex.hpp"
#include "sigpack.h"
void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[])
{
// Check input parameters
int M = mxGetM(prhs[0]);
int N = mxGetN(prhs[0]);
if(nrhs != 1) mexErrMsgTxt("One input required.");
if(nlhs != 1) mexErrMsgTxt("One output required.");
if( M == 1 && N >= 1 ) mexErrMsgTxt("Column vector required.");
if( mxIsComplex(prhs[0])) mexErrMsgTxt("Real input data required.");
// Convert each input argument to Armadillo
vec x = conv_to<vec>::from(armaGetPr(prhs[0],true));
//===================================================================
// START: Your stuff
int Nfft = 512;
vec y(x.size());
vec Pyy(Nfft);
int K = 17;
vec b(K);
// Create a FIR filter
b = sp::fir1(K,0.25);
fir_filt.set_coeffs(b);
// Filter the signal
y = fir_filt.filter(x);
// Calculate spectrum
Pyy = sp::pwelch(y,Nfft,Nfft/2);
// END: Your stuff
//===================================================================
// Convert outputs back to Matlab
plhs[0] = armaCreateMxMatrix(Pyy.size(), 1, mxDOUBLE_CLASS, mxREAL);
armaSetPr(plhs[0], conv_to<mat>::from(Pyy));
return;
}
FIR/MA filter class.
Definition filter.h:22
arma::Mat< T3 > filter(const arma::Mat< T1 > &in)
Filter function.
Definition filter.h:138
void set_coeffs(const arma::Mat< T2 > &_b)
Sets coefficients in FIR filter. The internal state and pointers are cleared.
Definition filter.h:70
arma_inline arma::vec fir1(const arma::uword M, const double f0)
FIR lowpass design function. FIR lowpassdesign using windows method (hamming window)....
Definition filter.h:637
arma::vec pwelch(const arma::Col< T1 > &x, const arma::uword Nfft=512, const arma::uword Noverl=256)
Power spectrum calculation using Welch's method.
Definition spectrum.h:166