getting ready for v100

This commit is contained in:
Mark Borgerding
2003-12-15 03:23:46 +00:00
parent 573536f48f
commit 6b23ebb5c1
13 changed files with 331 additions and 69 deletions

View File

@ -4,7 +4,6 @@ endif
ifeq "$(NUMFFTS)" ""
NUMFFTS=10000
endif
NROWS=30
ifeq "$(DATATYPE)" ""
DATATYPE=float
@ -32,7 +31,7 @@ CFLAGS=-Wall -O3
$(SELFTEST): ../kiss_fft.c $(SELFTESTSRC)
$(CC) -o $@ $(CFLAGS) -I.. $(TYPEFLAGS) -lm $+
$(TESTREAL): ../kiss_fft.c ../kiss_fftr.c test_real.c
$(TESTREAL): ../kiss_fft.c kiss_fftr.c test_real.c
$(CC) -o $@ $(CFLAGS) -I.. $(TYPEFLAGS) -lm $+
$(BENCHKISS): benchkiss.c ../kiss_fft.c pstats.c
@ -41,30 +40,24 @@ $(BENCHKISS): benchkiss.c ../kiss_fft.c pstats.c
$(BENCHFFTW): benchfftw.c pstats.c
@$(CC) -o $@ $(CFLAGS) -DDATATYPE$(DATATYPE) benchfftw.c pstats.c -lm -lfftw3f -lfftw3 -L /usr/local/lib/ || echo "FFTW not available for comparison"
POW2=256 512 1024 2048 4096 8192
POW3=243 729 2187
POW5=25 125 625 3125
mtime: all
@for n in $(POW2) $(POW3) $(POW5) ;do \
echo ============================;\
./$(BENCHKISS) -x $(NUMFFTS) -n $$n;\
[ -x ./$(BENCHFFTW) ] && ./$(BENCHFFTW) -x $(NUMFFTS) -n $$n || true ; \
done
test: all $(TESTREAL)
@echo "======SELF TEST"
test: all
@echo "======SELF TEST $(DATATYPE)"
@./$(SELFTEST)
@echo "======REAL FFT TEST"
@echo "======REAL FFT TEST $(DATATYPE)"
@./$(TESTREAL)
@echo "======TIMING TEST"
@echo "======TIMING TEST $(DATATYPE)"
@./$(BENCHKISS) -x $(NUMFFTS) -n $(NFFT)
@[ -x ./$(BENCHFFTW) ] && ./$(BENCHFFTW) -x $(NUMFFTS) -n $(NFFT) ||true
testall:
@(export DATATYPE=double && make test )
@(export DATATYPE=float && make test )
@(export DATATYPE=short && make test )
selftest.c:
./mk_test.py 10 12 14 > selftest.c
selftest_short.c:
./mk_test.py -s 10 12 14 > selftest_short.c
clean:
rm -f *~ bm_* st_* tr_* *.dat
rm -f *~ bm_* st_* tr_*

View File

@ -20,12 +20,6 @@ THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
#include "kiss_fft.h"
#ifndef NFFT
# define NFFT 1024
#endif
void fft_file(FILE * fin,FILE * fout,int nfft,int nrows,int isinverse,int useascii,int times)
{
int i;

140
tools/kiss_fftr.c Normal file
View File

@ -0,0 +1,140 @@
/*
Copyright (c) 2003, Mark Borgerding
All rights reserved.
Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:
* Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.
* Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution.
* Neither the author nor the names of any contributors may be used to endorse or promote products derived from this software without specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
*/
#include "_kiss_fft_guts.h"
typedef struct {
int minus3; /*magic to signify a 1-d real transform*/
kiss_fft_state * substate;
kiss_fft_cpx * tmpbuf;
kiss_fft_cpx * super_twiddles;
}kiss_fftr_state;
void * kiss_fftr_alloc(int nfft,int inverse_fft,void * mem,size_t * lenmem)
{
int i;
kiss_fftr_state *st = NULL;
size_t subsize, memneeded;
if (nfft & 1) {
/*fprintf(stderr,"Real FFT optimization must be even.\n"); */
return NULL;
}
nfft >>= 1;
kiss_fft_alloc (nfft, inverse_fft, NULL, &subsize);
memneeded = sizeof(kiss_fftr_state) + subsize + sizeof(kiss_fft_cpx) * ( nfft * 2);
if (lenmem == NULL) {
st = (kiss_fftr_state *) malloc (memneeded);
} else {
if (*lenmem >= memneeded)
st = (kiss_fftr_state *) mem;
*lenmem = memneeded;
}
if (!st)
return NULL;
st->minus3 = -3;
st->substate = (kiss_fft_state *) (st + 1); /*just beyond kiss_fftr_state struct */
st->tmpbuf = (kiss_fft_cpx *) (((char *) st->substate) + subsize);
st->super_twiddles = st->tmpbuf + nfft;
kiss_fft_alloc(nfft, inverse_fft, st->substate, &subsize);
for (i = 0; i < nfft; ++i) {
double phase =
-3.14159265358979323846264338327 * ((double) i / nfft + .5);
if (inverse_fft)
phase *= -1;
st->super_twiddles[i] = kf_cexp (phase);
}
return st;
}
void kiss_fftr(const void * cfg,const kiss_fft_scalar *timedata,kiss_fft_cpx *freqdata)
{
/* input buffer timedata is stored row-wise */
kiss_fftr_state *st = ( kiss_fftr_state *)cfg;
int k,N;
if ( st->minus3 != -3 || st->substate->inverse) {
fprintf(stderr,"kiss fft usage error: improper alloc\n");
exit(1);
}
N = st->substate->nfft;
/*perform the parallel fft of two real signals packed in real,imag*/
kiss_fft( st->substate , (const kiss_fft_cpx*)timedata, st->tmpbuf );
freqdata[0].r = st->tmpbuf[0].r + st->tmpbuf[0].i;
freqdata[0].i = 0;
C_FIXDIV(freqdata[0],2);
for (k=1;k <= N/2 ; ++k ) {
kiss_fft_cpx fpnk,fpk,f1k,f2k,tw;
fpk = st->tmpbuf[k];
fpnk.r = st->tmpbuf[N-k].r;
fpnk.i = -st->tmpbuf[N-k].i;
C_FIXDIV(fpk,2);
C_FIXDIV(fpnk,2);
C_ADD( f1k, fpk , fpnk );
C_SUB( f2k, fpk , fpnk );
C_MUL( tw , f2k , st->super_twiddles[k]);
C_ADD( freqdata[k] , f1k ,tw);
freqdata[k].r = (f1k.r + tw.r) / 2;
freqdata[k].i = (f1k.i + tw.i) / 2;
freqdata[N-k].r = (f1k.r - tw.r)/2;
freqdata[N-k].i = - (f1k.i - tw.i)/2;
}
freqdata[N].r = st->tmpbuf[0].r - st->tmpbuf[0].i;
freqdata[N].i = 0;
C_FIXDIV(freqdata[N],2);
}
void kiss_fftri(const void * cfg,const kiss_fft_cpx *freqdata,kiss_fft_scalar *timedata)
{
/* input buffer timedata is stored row-wise */
kiss_fftr_state *st = (kiss_fftr_state *) cfg;
int k, N;
if (st->minus3 != -3 || st->substate->inverse == 0) {
fprintf (stderr, "kiss fft usage error: improper alloc\n");
exit (1);
}
N = st->substate->nfft;
st->tmpbuf[0].r = freqdata[0].r + freqdata[N].r;
st->tmpbuf[0].i = freqdata[0].r - freqdata[N].r;
for (k = 1; k <= N / 2; ++k) {
kiss_fft_cpx fk, fnkc, fek, fok, tmpbuf;
fk = freqdata[k];
fnkc.r = freqdata[N - k].r;
fnkc.i = -freqdata[N - k].i;
C_ADD (fek, fk, fnkc);
C_SUB (tmpbuf, fk, fnkc);
C_MUL (fok, tmpbuf, st->super_twiddles[k]);
C_ADD (st->tmpbuf[k], fek, fok);
C_SUB (st->tmpbuf[N - k], fek, fok);
st->tmpbuf[N - k].i *= -1;
}
kiss_fft (st->substate, st->tmpbuf, (kiss_fft_cpx *) timedata);
}

12
tools/kiss_fftr.h Normal file
View File

@ -0,0 +1,12 @@
#ifndef KISS_FTR_H
#define KISS_FTR_H
#include "kiss_fft.h"
/* Real optimized version can save about 45% cpu time vs. complex fft of a real seq.
*/
void * kiss_fftr_alloc(int nfft,int inverse_fft,void * mem, size_t * lenmem);
void kiss_fftr(const void * cfg,const kiss_fft_scalar *timedata,kiss_fft_cpx *freqdata);
void kiss_fftri(const void * cfg,const kiss_fft_cpx *freqdata,kiss_fft_scalar *timedata);
#endif