123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557 |
- unit u_testfftunit;
- interface
- uses Math, Sysutils, u_ap, u_ftbase, u_fft;
- function TestFFT(Silent : Boolean):Boolean;
- function testfftunit_test_silent():Boolean;
- function testfftunit_test():Boolean;
- implementation
- procedure RefFFTC1D(var A : TComplex1DArray; N : Integer);forward;
- procedure RefFFTC1DInv(var A : TComplex1DArray; N : Integer);forward;
- procedure RefInternalCFFT(var a : TReal1DArray;
- nn : Integer;
- InverseFFT : Boolean);forward;
- procedure RefInternalRFFT(const a : TReal1DArray;
- nn : Integer;
- var F : TComplex1DArray);forward;
- (*************************************************************************
- Test
- *************************************************************************)
- function TestFFT(Silent : Boolean):Boolean;
- var
- N : Integer;
- I : Integer;
- K : Integer;
- A1 : TComplex1DArray;
- A2 : TComplex1DArray;
- A3 : TComplex1DArray;
- R1 : TReal1DArray;
- R2 : TReal1DArray;
- Buf : TReal1DArray;
- Plan : FTPlan;
- MaxN : Integer;
- BiDiErr : Double;
- BiDiRErr : Double;
- RefErr : Double;
- RefRErr : Double;
- REIntErr : Double;
- ErrTol : Double;
- RefErrors : Boolean;
- BiDiErrors : Boolean;
- RefRErrors : Boolean;
- BiDiRErrors : Boolean;
- REIntErrors : Boolean;
- WasErrors : Boolean;
- begin
- MaxN := 128;
- ErrTol := 100000*Power(MaxN, AP_Double(3)/2)*MachineEpsilon;
- BiDiErrors := False;
- RefErrors := False;
- BiDiRErrors := False;
- RefRErrors := False;
- REIntErrors := False;
- WasErrors := False;
-
- //
- // Test bi-directional error: norm(x-invFFT(FFT(x)))
- //
- BiDiErr := 0;
- BiDiRErr := 0;
- N:=1;
- while N<=MaxN do
- begin
-
- //
- // Complex FFT/invFFT
- //
- SetLength(A1, N);
- SetLength(A2, N);
- SetLength(A3, N);
- I:=0;
- while I<=N-1 do
- begin
- A1[I].X := 2*RandomReal-1;
- A1[I].Y := 2*RandomReal-1;
- A2[I] := A1[I];
- A3[I] := A1[I];
- Inc(I);
- end;
- FFTC1D(A2, N);
- FFTC1DInv(A2, N);
- FFTC1DInv(A3, N);
- FFTC1D(A3, N);
- I:=0;
- while I<=N-1 do
- begin
- BiDiErr := Max(BiDiErr, AbsComplex(C_Sub(A1[I],A2[I])));
- BiDiErr := Max(BiDiErr, AbsComplex(C_Sub(A1[I],A3[I])));
- Inc(I);
- end;
-
- //
- // Real
- //
- SetLength(R1, N);
- SetLength(R2, N);
- I:=0;
- while I<=N-1 do
- begin
- R1[I] := 2*RandomReal-1;
- R2[I] := R1[I];
- Inc(I);
- end;
- FFTR1D(R2, N, A1);
- APVMul(@R2[0], 0, N-1, 0);
- FFTR1DInv(A1, N, R2);
- I:=0;
- while I<=N-1 do
- begin
- BiDiRErr := Max(BiDiRErr, AbsComplex(C_Complex(R1[I]-R2[I])));
- Inc(I);
- end;
- Inc(N);
- end;
- BiDiErrors := BiDiErrors or AP_FP_Greater(BiDiErr,ErrTol);
- BiDiRErrors := BiDiRErrors or AP_FP_Greater(BiDiRErr,ErrTol);
-
- //
- // Test against reference O(N^2) implementation
- //
- RefErr := 0;
- RefRErr := 0;
- N:=1;
- while N<=MaxN do
- begin
-
- //
- // Complex FFT
- //
- SetLength(A1, N);
- SetLength(A2, N);
- I:=0;
- while I<=N-1 do
- begin
- A1[I].X := 2*RandomReal-1;
- A1[I].Y := 2*RandomReal-1;
- A2[I] := A1[I];
- Inc(I);
- end;
- FFTC1D(A1, N);
- RefFFTC1D(A2, N);
- I:=0;
- while I<=N-1 do
- begin
- RefErr := Max(RefErr, AbsComplex(C_Sub(A1[I],A2[I])));
- Inc(I);
- end;
-
- //
- // Complex inverse FFT
- //
- SetLength(A1, N);
- SetLength(A2, N);
- I:=0;
- while I<=N-1 do
- begin
- A1[I].X := 2*RandomReal-1;
- A1[I].Y := 2*RandomReal-1;
- A2[I] := A1[I];
- Inc(I);
- end;
- FFTC1DInv(A1, N);
- RefFFTC1DInv(A2, N);
- I:=0;
- while I<=N-1 do
- begin
- RefErr := Max(RefErr, AbsComplex(C_Sub(A1[I],A2[I])));
- Inc(I);
- end;
-
- //
- // Real forward/inverse FFT:
- // * calculate and check forward FFT
- // * use precalculated FFT to check backward FFT
- // fill unused parts of frequencies array with random numbers
- // to ensure that they are not really used
- //
- SetLength(R1, N);
- SetLength(R2, N);
- I:=0;
- while I<=N-1 do
- begin
- R1[I] := 2*RandomReal-1;
- R2[I] := R1[I];
- Inc(I);
- end;
- FFTR1D(R1, N, A1);
- RefInternalRFFT(R2, N, A2);
- I:=0;
- while I<=N-1 do
- begin
- RefRErr := Max(RefRErr, AbsComplex(C_Sub(A1[I],A2[I])));
- Inc(I);
- end;
- SetLength(A3, Floor(AP_Double(N)/2)+1);
- I:=0;
- while I<=Floor(AP_Double(N)/2) do
- begin
- A3[I] := A2[I];
- Inc(I);
- end;
- A3[0].Y := 2*RandomReal-1;
- if N mod 2=0 then
- begin
- A3[Floor(AP_Double(N)/2)].Y := 2*RandomReal-1;
- end;
- I:=0;
- while I<=N-1 do
- begin
- R1[I] := 0;
- Inc(I);
- end;
- FFTR1DInv(A3, N, R1);
- I:=0;
- while I<=N-1 do
- begin
- RefRErr := Max(RefRErr, AbsReal(R2[I]-R1[I]));
- Inc(I);
- end;
- Inc(N);
- end;
- RefErrors := RefErrors or AP_FP_Greater(RefErr,ErrTol);
- RefRErrors := RefRErrors or AP_FP_Greater(RefRErr,ErrTol);
-
- //
- // test internal real even FFT
- //
- REIntErr := 0;
- K:=1;
- while K<=MaxN div 2 do
- begin
- N := 2*K;
-
- //
- // Real forward FFT
- //
- SetLength(R1, N);
- SetLength(R2, N);
- I:=0;
- while I<=N-1 do
- begin
- R1[I] := 2*RandomReal-1;
- R2[I] := R1[I];
- Inc(I);
- end;
- FTBaseGenerateComplexFFTPlan(N div 2, Plan);
- SetLength(Buf, N);
- FFTR1DInternalEven(R1, N, Buf, Plan);
- RefInternalRFFT(R2, N, A2);
- REIntErr := Max(REIntErr, AbsReal(R1[0]-A2[0].X));
- REIntErr := Max(REIntErr, AbsReal(R1[1]-A2[N div 2].X));
- I:=1;
- while I<=N div 2-1 do
- begin
- REIntErr := Max(REIntErr, AbsReal(R1[2*I+0]-A2[I].X));
- REIntErr := Max(REIntErr, AbsReal(R1[2*I+1]-A2[I].Y));
- Inc(I);
- end;
-
- //
- // Real backward FFT
- //
- SetLength(R1, N);
- I:=0;
- while I<=N-1 do
- begin
- R1[I] := 2*RandomReal-1;
- Inc(I);
- end;
- SetLength(A2, Floor(AP_Double(N)/2)+1);
- A2[0] := C_Complex(R1[0]);
- I:=1;
- while I<=Floor(AP_Double(N)/2)-1 do
- begin
- A2[I].X := R1[2*I+0];
- A2[I].Y := R1[2*I+1];
- Inc(I);
- end;
- A2[Floor(AP_Double(N)/2)] := C_Complex(R1[1]);
- FTBaseGenerateComplexFFTPlan(N div 2, Plan);
- SetLength(Buf, N);
- FFTR1DInvInternalEven(R1, N, Buf, Plan);
- FFTR1DInv(A2, N, R2);
- I:=0;
- while I<=N-1 do
- begin
- REIntErr := Max(REIntErr, AbsReal(R1[I]-R2[I]));
- Inc(I);
- end;
- Inc(K);
- end;
- REIntErrors := REIntErrors or AP_FP_Greater(REIntErr,ErrTol);
-
- //
- // end
- //
- WasErrors := BiDiErrors or BiDiRErrors or RefErrors or RefRErrors or REIntErrors;
- if not Silent then
- begin
- Write(Format('TESTING FFT'#13#10'',[]));
- Write(Format('FINAL RESULT: ',[]));
- if WasErrors then
- begin
- Write(Format('FAILED'#13#10'',[]));
- end
- else
- begin
- Write(Format('OK'#13#10'',[]));
- end;
- Write(Format('* BI-DIRECTIONAL COMPLEX TEST: ',[]));
- if BiDiErrors then
- begin
- Write(Format('FAILED'#13#10'',[]));
- end
- else
- begin
- Write(Format('OK'#13#10'',[]));
- end;
- Write(Format('* AGAINST REFERENCE COMPLEX FFT: ',[]));
- if RefErrors then
- begin
- Write(Format('FAILED'#13#10'',[]));
- end
- else
- begin
- Write(Format('OK'#13#10'',[]));
- end;
- Write(Format('* BI-DIRECTIONAL REAL TEST: ',[]));
- if BiDiRErrors then
- begin
- Write(Format('FAILED'#13#10'',[]));
- end
- else
- begin
- Write(Format('OK'#13#10'',[]));
- end;
- Write(Format('* AGAINST REFERENCE REAL FFT: ',[]));
- if RefRErrors then
- begin
- Write(Format('FAILED'#13#10'',[]));
- end
- else
- begin
- Write(Format('OK'#13#10'',[]));
- end;
- Write(Format('* INTERNAL EVEN FFT: ',[]));
- if REIntErrors then
- begin
- Write(Format('FAILED'#13#10'',[]));
- end
- else
- begin
- Write(Format('OK'#13#10'',[]));
- end;
- if WasErrors then
- begin
- Write(Format('TEST FAILED'#13#10'',[]));
- end
- else
- begin
- Write(Format('TEST PASSED'#13#10'',[]));
- end;
- end;
- Result := not WasErrors;
- end;
- (*************************************************************************
- Reference FFT
- *************************************************************************)
- procedure RefFFTC1D(var A : TComplex1DArray; N : Integer);
- var
- Buf : TReal1DArray;
- I : Integer;
- begin
- Assert(N>0, 'FFTC1D: incorrect N!');
- SetLength(Buf, 2*N);
- I:=0;
- while I<=N-1 do
- begin
- Buf[2*I+0] := A[I].X;
- Buf[2*I+1] := A[I].Y;
- Inc(I);
- end;
- RefInternalCFFT(Buf, N, False);
- I:=0;
- while I<=N-1 do
- begin
- A[I].X := Buf[2*I+0];
- A[I].Y := Buf[2*I+1];
- Inc(I);
- end;
- end;
- (*************************************************************************
- Reference inverse FFT
- *************************************************************************)
- procedure RefFFTC1DInv(var A : TComplex1DArray; N : Integer);
- var
- Buf : TReal1DArray;
- I : Integer;
- begin
- Assert(N>0, 'FFTC1DInv: incorrect N!');
- SetLength(Buf, 2*N);
- I:=0;
- while I<=N-1 do
- begin
- Buf[2*I+0] := A[I].X;
- Buf[2*I+1] := A[I].Y;
- Inc(I);
- end;
- RefInternalCFFT(Buf, N, True);
- I:=0;
- while I<=N-1 do
- begin
- A[I].X := Buf[2*I+0];
- A[I].Y := Buf[2*I+1];
- Inc(I);
- end;
- end;
- (*************************************************************************
- Internal complex FFT stub.
- Uses straightforward formula with O(N^2) complexity.
- *************************************************************************)
- procedure RefInternalCFFT(var a : TReal1DArray;
- nn : Integer;
- InverseFFT : Boolean);
- var
- Tmp : TReal1DArray;
- I : Integer;
- J : Integer;
- K : Integer;
- Hre : Double;
- Him : Double;
- C : Double;
- S : Double;
- re : Double;
- im : Double;
- begin
- SetLength(Tmp, 2*nn-1+1);
- if not InverseFFT then
- begin
- i:=0;
- while i<=nn-1 do
- begin
- HRe := 0;
- Him := 0;
- K:=0;
- while K<=nn-1 do
- begin
- re := A[2*k];
- im := A[2*k+1];
- C := Cos(-2*Pi*k*i/nn);
- S := Sin(-2*Pi*k*i/nn);
- HRe := HRe+C*Re-S*im;
- HIm := HIm+C*im+S*re;
- Inc(K);
- end;
- Tmp[2*i] := Hre;
- Tmp[2*i+1] := Him;
- Inc(i);
- end;
- i:=0;
- while i<=2*nn-1 do
- begin
- A[I] := Tmp[I];
- Inc(i);
- end;
- end
- else
- begin
- k:=0;
- while k<=nn-1 do
- begin
- HRe := 0;
- Him := 0;
- i:=0;
- while i<=nn-1 do
- begin
- re := A[2*i];
- im := A[2*i+1];
- C := Cos(2*Pi*k*i/nn);
- S := Sin(2*Pi*k*i/nn);
- HRe := HRe+C*Re-S*im;
- HIm := HIm+C*im+S*re;
- Inc(i);
- end;
- Tmp[2*k] := Hre/nn;
- Tmp[2*k+1] := Him/nn;
- Inc(k);
- end;
- i:=0;
- while i<=2*nn-1 do
- begin
- A[I] := Tmp[I];
- Inc(i);
- end;
- end;
- end;
- (*************************************************************************
- Internal real FFT stub.
- Uses straightforward formula with O(N^2) complexity.
- *************************************************************************)
- procedure RefInternalRFFT(const a : TReal1DArray;
- nn : Integer;
- var F : TComplex1DArray);
- var
- Tmp : TReal1DArray;
- I : Integer;
- begin
- SetLength(Tmp, 2*nn-1+1);
- I:=0;
- while I<=nn-1 do
- begin
- Tmp[2*I] := A[I];
- Tmp[2*I+1] := 0;
- Inc(I);
- end;
- RefInternalCFFT(Tmp, nn, False);
- SetLength(F, nn);
- I:=0;
- while I<=nn-1 do
- begin
- F[I].X := Tmp[2*I+0];
- F[I].Y := Tmp[2*I+1];
- Inc(I);
- end;
- end;
- (*************************************************************************
- Silent unit test
- *************************************************************************)
- function testfftunit_test_silent():Boolean;
- begin
- Result := TestFFT(True);
- end;
- (*************************************************************************
- Unit test
- *************************************************************************)
- function testfftunit_test():Boolean;
- begin
- Result := TestFFT(False);
- end;
- end.
|