[3deeb5]: mex / comp_fftreal.c Maximize Restore History

Download this file

comp_fftreal.c    176 lines (127 with data), 3.8 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
#include "fftw3.h"
#include "stdlib.h"
#ifndef _LTFAT_MEX_FILE
#define _LTFAT_MEX_FILE
#define ISNARGINEQ 1
#define TYPEDEPARGS 0
#define SINGLEARGS
#define REALARGS
#define NOCOMPLEXFMTCHANGE
#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 "config.h"
static LTFAT_FFTW(plan)* LTFAT_NAME(p_old) = 0;
void LTFAT_NAME(fftrealAtExit)()
{
if(LTFAT_NAME(p_old)!=0)
{
LTFAT_FFTW(destroy_plan)(*LTFAT_NAME(p_old));
free(LTFAT_NAME(p_old));
}
}
// Calling convention:
// comp_fftreal(f);
void LTFAT_NAME(ltfatMexFnc)( int nlhs, mxArray *plhs[],int nrhs, const mxArray *prhs[] )
{
static int atExitRegistered = 0;
if(!atExitRegistered)
{
LTFAT_NAME(ltfatMexAtExit)(LTFAT_NAME(fftrealAtExit));
atExitRegistered = 1;
}
int L, W, L2;
LTFAT_REAL *f, *cout_r, *cout_i;
LTFAT_FFTW(iodim) dims[1], howmanydims[1];
LTFAT_FFTW(plan) p;
L = mxGetM(prhs[0]);
W = mxGetN(prhs[0]);
L2 = (L/2)+1;
// Get pointer to input.
f= (LTFAT_REAL*) mxGetPr(prhs[0]);
plhs[0] = ltfatCreateMatrix(L2, W, LTFAT_MX_CLASSID, mxCOMPLEX);
// Get pointer to output.
cout_r = (LTFAT_REAL*) mxGetPr(plhs[0]);
cout_i = (LTFAT_REAL*) mxGetPi(plhs[0]);
// Create plan. Copy data from f to cout.
dims[0].n = L;
dims[0].is = 1;
dims[0].os = 1;
howmanydims[0].n = W;
howmanydims[0].is = L;
howmanydims[0].os = L2;
// The calling prototype
//fftw_plan fftw_plan_guru_split_dft_r2c(
// int rank, const fftw_iodim *dims,
// int howmany_rank, const fftw_iodim *howmany_dims,
// double *in, double *ro, double *io,
// unsigned flags);
/* We are violating this here:
You must create the plan before initializing the input, because FFTW_MEASURE overwrites the in/out arrays.
(Technically, FFTW_ESTIMATE does not touch your arrays, but you should always create plans first just to be sure.)
*/
p = LTFAT_FFTW(plan_guru_split_dft_r2c)(1, dims,
1, howmanydims,
f, cout_r, cout_i, FFTW_ESTIMATE);
/*
FFTW documentation qote http://www.fftw.org/fftw3_doc/New_002darray-Execute-Functions.html#New_002darray-Execute-Functions:
...
creating a new plan is quick once one exists for a given size
...
so why not to store the old plan..
*/
LTFAT_NAME(fftrealAtExit)();
LTFAT_NAME(p_old) = malloc(sizeof(p));
memcpy(LTFAT_NAME(p_old),&p,sizeof(p));
// Real FFT.
LTFAT_FFTW(execute)(p);
// LTFAT_FFTW(destroy_plan)(p);
return;
}
#endif
/*
#include "mex.h"
#include "config.h"
#include "fftw3.h"
// Calling convention:
// comp_fftreal(f);
void mexFunction( int nlhs, mxArray *plhs[],
int nrhs, const mxArray *prhs[] )
{
int L, W, L2;
double *f, *cout_r, *cout_i;
fftw_iodim dims[1], howmanydims[1];
fftw_plan p;
L = mxGetM(prhs[0]);
W = mxGetN(prhs[0]);
L2 = (L/2)+1;
// Get pointer to input.
f=mxGetPr(prhs[0]);
plhs[0] = mxCreateDoubleMatrix(L2, W, mxCOMPLEX);
// Get pointer to output.
cout_r = mxGetPr(plhs[0]);
cout_i = mxGetPi(plhs[0]);
// Create plan. Copy data from f to cout.
dims[0].n = L;
dims[0].is = 1;
dims[0].os = 1;
howmanydims[0].n = W;
howmanydims[0].is = L;
howmanydims[0].os = L2;
// The calling prototype
//fftw_plan fftw_plan_guru_split_dft_r2c(
// int rank, const fftw_iodim *dims,
// int howmany_rank, const fftw_iodim *howmany_dims,
// double *in, double *ro, double *io,
// unsigned flags);
p = fftw_plan_guru_split_dft_r2c(1, dims,
1, howmanydims,
f, cout_r, cout_i, FFTW_OPTITYPE);
// Real FFT.
fftw_execute(p);
fftw_destroy_plan(p);
return;
}
*/