/* Copyright (C) 2004 Caltech Written by Lexing Ying */ #include "fdct_usfft.hpp" #include "fdct_usfft_inline.hpp" FDCT_USFFT_NS_BEGIN_NAMESPACE //------------------------------------ int fdct_usfft_fftL(int N1, int N2, CpxNumMat& x, CpxOffMat& O); int fdct_usfft_sepscale(int N1, int N2, int nbscales, int ac, CpxOffMat& O, vector& Xhghs); int fdct_usfft_sepangle(double XL1, double XL2, int nbangles, CpxOffMat& Xhgh, vector& msc); int fdct_usfft_ifftS(vector& msc, vector& csc); int fdct_usfft_wavelet(CpxOffMat& Xhgh, vector& csc); int fdct_usfft_1dinterp(DblNumMat& off, DblNumMat& wgt, CpxOffVec& val, CpxNumMat& res, map& f1map, map& b1map); //------------------------------------ int fdct_usfft( int N1, int N2, int nbscales, int nbangles_coarse, int ac, CpxNumMat& x, vector< vector >& c) { assert(N1==x.m() && N2==x.n()); int F1 = N1/2; int F2 = N2/2; //1. fftL CpxOffMat O(N1, N2); fdct_usfft_fftL(N1, N2, x, O); //2. seperate scale vector Xhghs(nbscales); fdct_usfft_sepscale(N1, N2, nbscales, ac, O, Xhghs); //3. work on each scale if(ac==1) { //nbangles vector nbangles(nbscales); nbangles[0] = 1; for(int sc=1; sc0; sc--) { vector msc(nbangles[sc]); fdct_usfft_sepangle(XL1, XL2, nbangles[sc], Xhghs[sc], msc); fdct_usfft_ifftS(msc, c[sc]); XL1 = XL1/2; XL2 = XL2/2; } //coarest level fdct_usfft_wavelet(Xhghs[0], c[0]); } else { //nbangles vector nbangles(nbscales); nbangles[0] = 1; for(int sc=1; sc0; sc--) { vector msc(nbangles[sc]); fdct_usfft_sepangle(XL1, XL2, nbangles[sc], Xhghs[sc], msc); fdct_usfft_ifftS(msc, c[sc]); XL1 = XL1/2; XL2 = XL2/2; } //coarest level fdct_usfft_wavelet(Xhghs[0], c[0]); } return 0; } //-------------------- int fdct_usfft_fftL(int N1, int N2, CpxNumMat& x, CpxOffMat& O) { //move zero freq at front int F1 = N1/2; int F2 = N2/2; //1. ifftshift CpxNumMat T(x); //2. do fft and scale fftwnd_plan p = fftw2d_create_plan(N2, N1, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_IN_PLACE); fftwnd_one(p, (fftw_complex*)T.data(), NULL); fftwnd_destroy_plan(p); double sqrtprod = sqrt(double(N1*N2)); for(int i=0; i& Xhghs) { //----------------------------------------------------- //unfold if necessary CpxOffMat X; if(ac==1) { //-------------------------- double XL1 = 4.0*N1/3.0; double XL2 = 4.0*N2/3.0; //range int XS1, XS2; int XF1, XF2; double XR1, XR2; fdct_usfft_rangecompute(XL1, XL2, XS1, XS2, XF1, XF2, XR1, XR2); IntOffVec t1(XS1); for(int i=-XF1; i<-XF1+XS1; i++) if( i<-N1/2) t1(i) = i+int(N1); else if(i>(N1-1)/2) t1(i) = i-int(N1); else t1(i) = i; IntOffVec t2(XS2); for(int i=-XF2; i<-XF2+XS2; i++) if( i<-N2/2) t2(i) = i+int(N2); else if(i>(N2-1)/2) t2(i) = i-int(N2); else t2(i) = i; X.resize(XS1, XS2); for(int i=-XF1; i<-XF1+XS1; i++) for(int j=-XF2; j<-XF2+XS2; j++) X(i,j) = O(t1(i), t2(j)); DblOffMat lowpass(XS1,XS2); fdct_usfft_lowpasscompute(XL1, XL2, lowpass); //compute the low pass filter for(int i=-XF1; i<-XF1+XS1; i++) for(int j=-XF2; j<-XF2+XS2; j++) X(i,j) *= lowpass(i,j); } else { //-------------------------- X = O; } //----------------------------------------------------- //seperate Xhghs.resize(nbscales); double XL1 = 4.0*N1/3.0; double XL2 = 4.0*N2/3.0; //range //int XS1, XS2; int XF1, XF2; double XR1, XR2; fdct_usfft_rangecompute(XL1, XL2, XS1, XS2, XF1, XF2, XR1, XR2); for(int sc=nbscales-1; sc>0; sc--) { double XL1n = XL1/2; double XL2n = XL2/2; int XS1n, XS2n; int XF1n, XF2n; double XR1n, XR2n; fdct_usfft_rangecompute(XL1n, XL2n, XS1n, XS2n, XF1n, XF2n, XR1n, XR2n); //get filters DblOffMat lowpass(XS1n, XS2n); fdct_usfft_lowpasscompute(XL1n, XL2n, lowpass); DblOffMat hghpass(XS1n, XS2n); for(int i=-XF1n; i<-XF1n+XS1n; i++) for(int j=-XF2n; j<-XF2n+XS2n; j++) hghpass(i,j) = sqrt(1-lowpass(i,j)*lowpass(i,j)); //get Xhgh CpxOffMat Xhgh(X); for(int i=-XF1n; i<-XF1n+XS1n; i++) for(int j=-XF2n; j<-XF2n+XS2n; j++) { //Xhgh(i,j).re *= hghpass(i,j); Xhgh(i,j).im *= hghpass(i,j); Xhgh(i,j) *= hghpass(i,j); } CpxOffMat Xlow(XS1n, XS2n); for(int i=-XF1n; i<-XF1n+XS1n; i++) for(int j=-XF2n; j<-XF2n+XS2n; j++) { //Xlow(i,j).re = X(i,j).re * lowpass(i,j); Xlow(i,j).im = X(i,j).im * lowpass(i,j); Xlow(i,j) = X(i,j) * lowpass(i,j); } //set into vector Xhghs[sc] = Xhgh; X = Xlow; XL1 = XL1/2; XL2 = XL2/2; } Xhghs[0] = X; return 0; } //-------------------- int fdct_usfft_sepangle(double XL1, double XL2, int nbangles, CpxOffMat& Xhgh, vector& msc) { //int XS1, XS2; int XF1, XF2; double XR1, XR2; fdct_usfft_rangecompute(XL1, XL2, XS1, XS2, XF1, XF2, XR1, XR2); map f1map; //forward 1d map, IN_PLACE map b1map; //backward 1d map, IN_PLACE //allocate msc msc.resize(nbangles); int nbquadrants = 4; int nd = nbangles / 4; int wcnt = 0; CpxOffMat Xhghb(Xhgh); double XL1b = XL1; double XL2b = XL2; int qvec[] = {2,1,0,3}; for(int qi=0; qi tsc(nd); for(int w=0; w=-xh+xn) tmpx-=xn; for(int w=0; w=0; w--) { fdct_usfft_rotate_backward(q, tsc[w], msc[wcnt]); wcnt++; } } //for loop for quadrant Xhgh = Xhghb; XL1 = XL1b; XL2 = XL2b; assert(wcnt==nbangles); for(map::iterator mit=f1map.begin(); mit!=f1map.end(); mit++) { fftw_plan p = (*mit).second; fftw_destroy_plan(p); } for(map::iterator mit=b1map.begin(); mit!=b1map.end(); mit++) { fftw_plan p = (*mit).second; fftw_destroy_plan(p); } return 0; } //-------------------- int fdct_usfft_ifftS(vector& msc, vector& csc) { typedef pair intpair; map planmap; //do work csc.resize(msc.size()); for(int w=0; w::iterator mit=planmap.find( intpair(xn,yn) ); fftwnd_plan p = NULL; if(mit!=planmap.end()) { p = (*mit).second; } else { p = fftw2d_create_plan(yn, xn, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_IN_PLACE); planmap[ intpair(xn,yn) ] = p; } fftwnd_one(p, (fftw_complex*)tpdata.data(), NULL); double sqrtprod = sqrt(double(xn*yn)); for(int i=0; i::iterator mit=planmap.begin(); mit!=planmap.end(); mit++) { fftwnd_plan p = (*mit).second; fftwnd_destroy_plan(p); } return 0; } //-------------------- int fdct_usfft_wavelet(CpxOffMat& Xhgh, vector& csc) { int N1 = Xhgh.m(); int N2 = Xhgh.n(); int F1 = -Xhgh.s(); int F2 = -Xhgh.t(); CpxNumMat T(N1, N2); fdct_usfft_ifftshift(Xhgh, T); fftwnd_plan p = fftw2d_create_plan(N2, N1, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_IN_PLACE); fftwnd_one(p, (fftw_complex*)T.data(), NULL); fftwnd_destroy_plan(p); double sqrtprod = sqrt(double(N1*N2)); for(int i=0; i& f1map, map& b1map) { if(off.n()<=64) { //SIMPLE INTERPOLATION //-------------------------------------------- int N = val.m(); int F = -val.s(); fftw_plan fp = NULL; map::iterator fit = f1map.find(N); if(fit!=f1map.end()) { fp = (*fit).second; } else { fp = fftw_create_plan(N, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_IN_PLACE); f1map[N] = fp; } fftw_plan bp = NULL; map::iterator bit = b1map.find(N); if(bit!=b1map.end()) { bp = (*bit).second; } else { bp = fftw_create_plan(N, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_IN_PLACE); b1map[N] = bp; } //CpxOffVec val //right value CpxOffVec feq(N); //right freq CpxOffVec vvv(N); //new value CpxOffVec fff(N); //new freq CpxNumVec tmp(N); //1. fft_mid val fdct_usfft_ifftshift(val, tmp); fftw_one(fp, (fftw_complex*)tmp.data(), NULL); //fft double scale = sqrt(double(N)); for(int k=0; kF) ? k-N : k; res(k,w) = vvv(rk) * wgt(k,w); } } } else { //USFFT INTERPOLATION //-------------------------------------------- int L = 4; int D = 16; int N = val.m(); int F = -val.s(); fftw_plan fp = NULL; map::iterator fit=f1map.find(N); if(fit!=f1map.end()) { fp = (*fit).second; } else { fp = fftw_create_plan(N, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_IN_PLACE); f1map[N] = fp; } fftw_plan bp = NULL; map::iterator bit = b1map.find(D*N); if(bit!=b1map.end()) { bp = (*bit).second; } else { bp = fftw_create_plan(D*N, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_IN_PLACE); b1map[D*N] = bp; } CpxNumVec tmp(N); CpxNumVec exttmp(D*N); CpxOffVec feq(N); vector extdat(L); //1. get freq fdct_usfft_ifftshift(val, tmp); fftw_one(fp, (fftw_complex*)tmp.data(), NULL); //here tmp are the real frequencies double scale = double(N); for(int k=0; k=M_PI) cof-= 2.0*M_PI; //collapse into [-pi,pi) int ind = int(floor(cof/step)); ind = min(max(ind, -D*N/2), D*N/2-1); double dta = cof-ind*step; double pow = 1; cpx sum(0.0,0.0);// for(int l=0; l