/* 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_adjfftL(int N1, int N2, CpxOffMat& O, CpxNumMat& x); int fdct_usfft_adjsepscale(int N1, int N2, int nbscales, int ac, vector& Xhghs, CpxOffMat& O); int fdct_usfft_adjsepangle(double XL1, double XL2, int nbangles, vector& msc, CpxOffMat& Xhgh); int fdct_usfft_adjifftS(vector& csc, vector& msc); int fdct_usfft_adjwavelet(vector& csc, CpxOffMat& Xhgh); int fdct_usfft_adj1dinterp(DblNumMat& off, DblNumMat& wgt, CpxNumMat& res, CpxOffVec& val, map& f1map, map& b1map); //------------------------------------ int afdct_usfft(int N1, int N2, int nbscales, int nbangles_coarse, int ac, vector< vector >& c, CpxNumMat& x) { assert(nbscales==c.size() && nbangles_coarse==c[1].size()); //int F1 = N1/2; int F2 = N2/2; //1. vector Xhghs; Xhghs.resize(nbscales); if(ac==1) { vector nbangles(nbscales); nbangles[0] = 1; for(int sc=1; sc0; sc--) { vector msc(nbangles[sc]); fdct_usfft_adjifftS(c[sc], msc); fdct_usfft_adjsepangle(XL1, XL2, nbangles[sc], msc, Xhghs[sc]); XL1 = XL1/2; XL2 = XL2/2; } fdct_usfft_adjwavelet(c[0], Xhghs[0]); } else { vector nbangles(nbscales); nbangles[0] = 1; for(int sc=1; sc0; sc--) { vector msc(nbangles[sc]); fdct_usfft_adjifftS(c[sc], msc); fdct_usfft_adjsepangle(XL1, XL2, nbangles[sc], msc, Xhghs[sc]); XL1 = XL1/2; XL2 = XL2/2; } fdct_usfft_adjwavelet(c[0], Xhghs[0]); } //2. CpxOffMat O; fdct_usfft_adjsepscale(N1, N2, nbscales, ac, Xhghs, O); //3. fdct_usfft_adjfftL(N1, N2, O, x); return 0; } //------------------------------------ int fdct_usfft_adjfftL(int N1, int N2, CpxOffMat& O, CpxNumMat& x) { int F1 = N1/2; int F2 = N2/2; //ifftshift CpxNumMat T(N1, N2); fdct_usfft_ifftshift(O, T); //ifft 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& Xhghs, CpxOffMat& O) { int F1 = N1/2; int F2 = N2/2; 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); X.resize(XS1, XS2); } else { X.resize(N1, N2); } double XL1 = 4.0*N1/3.0; double XL2 = 4.0*N2/3.0; 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); 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)); for(int i=-XF1n; i<-XF1n+XS1n; i++) for(int j=-XF2n; j<-XF2n+XS2n; j++) Xhghs[sc](i,j) *= hghpass(i,j); for(int i=-XF1n; i<-XF1n+XS1n; i++) for(int j=-XF2n; j<-XF2n+XS2n; j++) Xhghs[sc-1](i,j) *= lowpass(i,j); CpxOffMat& G = Xhghs[sc]; for(int i=G.s(); 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; for(int i=-XF1; i<-XF1+XS1; i++) for(int j=-XF2; j<-XF2+XS2; j++) O(t1(i), t2(j)) += X(i,j); } else { O = X; } return 0; } //------------------------------------ int fdct_usfft_adjsepangle(double XL1, double XL2, int nbangles, vector& msc, CpxOffMat& Xhgh) { 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 space Xhgh.resize(XS1, XS2); clear(Xhgh); assert(msc.size()==nbangles); int nd = nbangles / 4; int nbquadrants = 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); //rotate data forward into msc for(int w=nd-1; w>=0; w--) { fdct_usfft_rotate_forward(q, msc[wcnt], tsc[w]); wcnt++; } //the adjoint of sample using USFFT for each line for(int xcur=xf; xcur=-xh+xn) tmpx-=xn; CpxNumMat res(yn, nd); for(int w=0; w::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_adjifftS(vector& csc, vector& msc) { typedef pair intpair; map planmap; //do work msc.resize(csc.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_FORWARD, 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_adjwavelet(vector& csc, CpxOffMat& Xhgh) { assert(csc.size()==1); CpxNumMat& C = csc[0]; int N1 = C.m(); int N2 = C.n(); CpxNumMat T(C); //ifftshift //CpxNumMat T(N1, N2); fdct_usfft_ifftshift(N1, N2, F1, F2, csc[0], T); //fft 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& 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; } map::iterator bit = b1map.find(N); fftw_plan bp = NULL; if(bit!=b1map.end()) { bp = (*bit).second; } else { bp = fftw_create_plan(N, FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_IN_PLACE); b1map[N] = bp; } CpxOffVec feq(N); CpxOffVec vvv(N); CpxOffVec fff(N); CpxNumVec tmp(N); //1. for each wedge, from res, fft, multiply, add back int nbwedges = off.n(); CpxOffVec sum(N); //sum, equal to zero for(int w=0; wF)? k-N : k; vvv(rk) = res(k,w) * wgt(k,w); } //new value fdct_usfft_ifftshift(vvv, tmp); fftw_one(fp, (fftw_complex*)tmp.data(), NULL); double scale = sqrt(double(N)); for(int k=0; k::iterator fit=f1map.find(D*N); if(fit!=f1map.end()) { fp = (*fit).second; } else { fp = fftw_create_plan(D*N, FFTW_FORWARD, FFTW_ESTIMATE | FFTW_IN_PLACE); f1map[D*N] = fp; } fftw_plan bp = NULL; //length N 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; } CpxNumVec tmp(N); CpxNumVec exttmp(D*N); CpxOffVec feq(N); clear(feq); CpxOffVec pok(N); setvalue(pok, cpx(1,0)); //power of k, unit vector vector extdat(L); //CpxOffVec extdat[4]; for(int l=0; l=M_PI) cof-= 2*M_PI; //collapse into [-pi,pi) int ind = int(floor(cof/step)); ind = min(max(ind, -D*N/2), D*N/2-1); //assert(ind>=-D*N/2 && ind