LTFAT_EXTERN void LTFAT_NAME(gabreassign)(const LTFAT_REAL *s, const LTFAT_REAL *tgrad, const LTFAT_REAL *fgrad, const ltfatInt L, const ltfatInt W, const ltfatInt a, const ltfatInt M, LTFAT_REAL *sr) { ltfatInt ii, posi, posj; const ltfatInt N=L/a; const ltfatInt b=L/M; ltfatInt *timepos = ltfat_malloc(N*sizeof*timepos); ltfatInt *freqpos = ltfat_malloc(M*sizeof*freqpos); fftindex(N,timepos); fftindex(M,freqpos); /* Zero the output array. */ memset(sr,0,M*N*W*sizeof*sr); for (ltfatInt w=0; w<W; w++) { for (ii=0; ii<M; ii++) { for (ltfatInt jj=0; jj<N; jj++) { /* Do a 'round' followed by a 'mod'. 'round' is not * present in all libraries, so use trunc(x+.5) instead */ /*posi=positiverem((ltfatInt)trunc(tgrad[ii+jj*M]/b+freqpos[ii]+.5),M); posj=positiverem((ltfatInt)trunc(fgrad[ii+jj*M]/a+timepos[jj]+.5),N);*/ posi=positiverem(ltfat_round(tgrad[ii+jj*M]/b+freqpos[ii]),M); posj=positiverem(ltfat_round(fgrad[ii+jj*M]/a+timepos[jj]),N); sr[posi+posj*M]+=s[ii+jj*M]; } } } LTFAT_SAFEFREEALL(freqpos,timepos); }
LTFAT_EXTERN void LTFAT_NAME(idgt_fb)(const LTFAT_COMPLEX *cin, const LTFAT_COMPLEX *g, const int L, const int gl, const int W, const int a, const int M, LTFAT_COMPLEX *f) { /* --------- initial declarations -------------- */ const int N=L/a; int ep, sp; /* This is a floor operation. */ const int glh=gl/2; /* This is a ceil operation. */ const int glh_d_a=(int)ceil((glh*1.0)/(a)); LTFAT_COMPLEX *fw; LTFAT_COMPLEX *cbuf = (LTFAT_COMPLEX*)ltfat_malloc(M*sizeof(LTFAT_COMPLEX)); /* Create plan. In-place. */ LTFAT_FFTW(plan) p_small = LTFAT_FFTW(plan_dft_1d)(M, cbuf, cbuf, FFTW_BACKWARD, FFTW_MEASURE); /* % The fftshift actually makes some things easier. */ LTFAT_COMPLEX *gw = (LTFAT_COMPLEX*)ltfat_malloc(gl*sizeof(LTFAT_COMPLEX)); for (int l=0;l<glh;l++) { gw[l][0] = g[l+(gl-glh)][0]; gw[l][1] = g[l+(gl-glh)][1]; } for (int l=glh;l<gl;l++) { gw[l][0] = g[l-glh][0]; gw[l][1] = g[l-glh][1]; } LTFAT_COMPLEX *ff = (LTFAT_COMPLEX*)ltfat_malloc(gl*sizeof(LTFAT_COMPLEX)); for (int w=0; w<W; w++) { fw=f+w*L; for (int l=0;l<L;l++) { fw[l][0]=0.0; fw[l][1]=0.0; } /* ----- Handle the first boundary using periodic boundary conditions. --- */ for (int n=0; n<glh_d_a; n++) { THE_SUM; sp=positiverem(n*a-glh,L); ep=positiverem(n*a-glh+gl-1,L); /* % Add the ff vector to f at position sp. */ for (int ii=0;ii<L-sp;ii++) { fw[sp+ii][0]+=ff[ii][0]; fw[sp+ii][1]+=ff[ii][1]; } for (int ii=0; ii<ep+1;ii++) { fw[ii][0] +=ff[L-sp+ii][0]; fw[ii][1] +=ff[L-sp+ii][1]; } } /* ----- Handle the middle case. --------------------- */ for (int n=glh_d_a; n<(L-(gl+1)/2)/a+1; n++) { THE_SUM; sp=positiverem(n*a-glh,L); ep=positiverem(n*a-glh+gl-1,L); /* Add the ff vector to f at position sp. */ for (int ii=0;ii<ep-sp+1;ii++) { fw[ii+sp][0] += ff[ii][0]; fw[ii+sp][1] += ff[ii][1]; } } /* Handle the last boundary using periodic boundary conditions. */ for (int n=(L-(gl+1)/2)/a+1; n<N; n++) { THE_SUM; sp=positiverem(n*a-glh,L); ep=positiverem(n*a-glh+gl-1,L); /* Add the ff vector to f at position sp. */ for (int ii=0;ii<L-sp;ii++) { fw[sp+ii][0]+=ff[ii][0]; fw[sp+ii][1]+=ff[ii][1]; } for (int ii=0; ii<ep+1;ii++) { fw[ii][0] +=ff[L-sp+ii][0]; fw[ii][1] +=ff[L-sp+ii][1]; } } } ltfat_free(cbuf); ltfat_free(ff); ltfat_free(gw); LTFAT_FFTW(destroy_plan)(p_small); }
/* wfac for real valued input. */ LTFAT_EXTERN void LTFAT_NAME(wfac_r)(const LTFAT_REAL *g, const int L, const int R, const int a, const int M, LTFAT_COMPLEX *gf) { int h_a, h_m; LTFAT_REAL *sbuf, *gfp; int s; int rem, negrem; LTFAT_FFTW(plan) p_before; const int b=L/M; const int c=gcd(a, M,&h_a, &h_m); const int p=a/c; const int q=M/c; const int d=b/p; const double sqrtM=sqrt(M); sbuf = (LTFAT_REAL*)ltfat_malloc(2*d*sizeof(LTFAT_REAL)); /* Create plan. In-place. */ p_before = LTFAT_FFTW(plan_dft_1d)(d, (LTFAT_COMPLEX*)sbuf, (LTFAT_COMPLEX*)sbuf, FFTW_FORWARD, FFTW_MEASURE); const int ld3=c*p*q*R; gfp=(LTFAT_REAL*)gf; for (int r=0;r<c;r++) { for (int w=0;w<R;w++) { for (int l=0;l<q;l++) { for (int k=0;k<p;k++) { negrem = positiverem(k*M-l*a,L); for (s=0;s<d;s++) { rem = (negrem+s*p*M)%L; sbuf[2*s] = sqrtM*g[r+rem+L*w]; sbuf[2*s+1] = 0.0; } LTFAT_FFTW(execute)(p_before); for (s=0;s<2*d;s+=2) { gfp[s*ld3] = sbuf[s]; gfp[s*ld3+1]= sbuf[s+1]; } gfp+=2; } } } } ltfat_free(sbuf); }
void dgt_fac(ltfat_complex *f, ltfat_complex *gf, const int L, const int W, const int R, const int a, const int M, ltfat_complex *cout, int dotime) { /* --------- initial declarations -------------- */ int b, N, c, d, p, q, h_a, h_m; double *gbase, *fbase, *cbase; int l, k, r, s, u, w, rw, nm, mm, km; int ld2a, ld1b, ld3b; int ld4c, ld5c; int rem; fftw_plan p_before, p_after, p_veryend; double *ff, *cf, *ffp, *sbuf, *fp, *cfp; double scalconst; double st0, st1, st6, st7; /* ----------- calculation of parameters and plans -------- */ if (dotime) { st0=ltfat_time(); } b=L/M; N=L/a; c=gcd(a, M,&h_a, &h_m); p=a/c; q=M/c; d=b/p; h_a=-h_a; /* Scaling constant needed because of FFTWs normalization. */ scalconst=1.0/((double)d*sqrt((double)M)); ff = (double*)ltfat_malloc(2*d*p*q*W*sizeof(double)); cf = (double*)ltfat_malloc(2*d*q*q*W*R*sizeof(double)); sbuf = (double*)ltfat_malloc(2*d*sizeof(double)); /* Create plans. In-place. */ p_before = fftw_plan_dft_1d(d, (ltfat_complex*)sbuf, (ltfat_complex*)sbuf, FFTW_FORWARD, FFTW_MEASURE); p_after = fftw_plan_dft_1d(d, (ltfat_complex*)sbuf, (ltfat_complex*)sbuf, FFTW_BACKWARD, FFTW_MEASURE); /* Create plan. In-place. */ p_veryend = fftw_plan_many_dft(1, &M, N*R*W, cout, NULL, 1, M, cout, NULL, 1, M, FFTW_FORWARD, FFTW_OPTITYPE); if (dotime) { st1=ltfat_time(); printf("DGT_FAC_7: Planning phase %f\n",st1-st0); } /* Leading dimensions of the 4dim array. */ ld2a=2*p*q*W; /* Leading dimensions of cf */ ld1b=q*R; ld3b=2*q*R*q*W; /* --------- main loop begins here ------------------- */ for (r=0;r<c;r++) { /* ---------- compute signal factorization ----------- */ ffp=ff; fp=(double*)f+2*r; if (p==1) /* Integer oversampling case */ { for (w=0;w<W;w++) { for (l=0;l<q;l++) { for (s=0;s<d;s++) { rem = 2*((s*M+l*a)%L); sbuf[2*s] = fp[rem]; sbuf[2*s+1] = fp[rem+1]; } fftw_execute(p_before); for (s=0;s<d;s++) { ffp[s*ld2a] = sbuf[2*s]*scalconst; ffp[s*ld2a+1] = sbuf[2*s+1]*scalconst; } ffp+=2; } fp+=2*L; } fp-=2*L*W; } else { /* rational sampling case */ for (w=0;w<W;w++) { for (l=0;l<q;l++) { for (k=0;k<p;k++) { for (s=0;s<d;s++) { rem = 2*positiverem(k*M+s*p*M-l*h_a*a, L); sbuf[2*s] = fp[rem]; sbuf[2*s+1] = fp[rem+1]; } fftw_execute(p_before); for (s=0;s<d;s++) { ffp[s*ld2a] = sbuf[2*s]*scalconst; ffp[s*ld2a+1] = sbuf[2*s+1]*scalconst; } ffp+=2; } } fp+=2*L; } fp-=2*L*W; } /* ----------- compute matrix multiplication ----------- */ /* Do the matmul */ if (p==1) { /* Integer oversampling case */ /* Rational oversampling case */ for (s=0;s<d;s++) { gbase=(double*)gf+2*(r+s*c)*q*R; fbase=ff+2*s*q*W; cbase=cf+2*s*q*q*W*R; for (nm=0;nm<q*W;nm++) { for (mm=0;mm<q*R;mm++) { cbase[0]=gbase[0]*fbase[0]+gbase[1]*fbase[1]; cbase[1]=gbase[0]*fbase[1]-gbase[1]*fbase[0]; gbase+=2; cbase+=2; } gbase-=2*q*R; fbase+=2; } cbase-=2*q*R*q*W; } } else { /* Rational oversampling case */ for (s=0;s<d;s++) { gbase=(double*)gf+2*(r+s*c)*p*q*R; fbase=ff+2*s*p*q*W; cbase=cf+2*s*q*q*W*R; for (nm=0;nm<q*W;nm++) { for (mm=0;mm<q*R;mm++) { cbase[0]=0.0; cbase[1]=0.0; for (km=0;km<p;km++) { cbase[0]+=gbase[0]*fbase[0]+gbase[1]*fbase[1]; cbase[1]+=gbase[0]*fbase[1]-gbase[1]*fbase[0]; gbase+=2; fbase+=2; } fbase-=2*p; cbase+=2; } gbase-=2*q*R*p; fbase+=2*p; } cbase-=2*q*R*q*W; fbase-=2*p*q*W; } } /* ------- compute inverse coefficient factorization ------- */ cfp=cf; ld4c=M*N; ld5c=M*N*R; /* Cover both integer and rational sampling case */ for (w=0;w<W;w++) { /* Complete inverse fac of coefficients */ for (l=0;l<q;l++) { for (rw=0;rw<R;rw++) { for (u=0;u<q;u++) { for (s=0;s<d;s++) { sbuf[2*s] = cfp[s*ld3b]; sbuf[2*s+1] = cfp[s*ld3b+1]; } cfp+=2; /* Do inverse fft of length d */ fftw_execute(p_after); for (s=0;s<d;s++) { rem= r+l*c+positiverem(u+s*q-l*h_a,N)*M+rw*ld4c+w*ld5c; cout[rem][0]=sbuf[2*s]; cout[rem][1]=sbuf[2*s+1]; } } } } } /* ----------- Main loop ends here ------------------------ */ } if (dotime) { st6=ltfat_time(); printf("DGT_FAC_7: Main loop done %f\n",st6-st1); } /* FFT to modulate the coefficients. */ fftw_execute(p_veryend); if (dotime) { st7=ltfat_time(); printf("DGT_FAC_7: Final FFT %f\n",st7-st6); printf("DGT_FAC_7: Total time %f\n",st7-st0); } /* ----------- Clean up ----------------- */ fftw_destroy_plan(p_before); fftw_destroy_plan(p_after); fftw_destroy_plan(p_veryend); ltfat_free(sbuf); ltfat_free(ff); ltfat_free(cf); }
/* wfac for real valued input. Produces only half the output coefficients of wfac_r */ LTFAT_EXTERN void LTFAT_NAME(wfacreal)(const LTFAT_REAL *g, const int L, const int R, const int a, const int M, LTFAT_COMPLEX *gf) { int h_a, h_m; LTFAT_REAL *gfp; int s; int rem, negrem; LTFAT_FFTW(plan) p_before; const int b=L/M; const int c=gcd(a, M,&h_a, &h_m); const int p=a/c; const int q=M/c; const int d=b/p; /* This is a floor operation. */ const int d2= d/2+1; const double sqrtM=sqrt(M); LTFAT_REAL *sbuf = (LTFAT_REAL*)ltfat_malloc(d*sizeof(LTFAT_REAL)); LTFAT_COMPLEX *cbuf = (LTFAT_COMPLEX*)ltfat_malloc(d2*sizeof(LTFAT_COMPLEX)); /* Create plan. In-place. */ p_before = LTFAT_FFTW(plan_dft_r2c_1d)(d, sbuf, cbuf, FFTW_MEASURE); const int ld3=2*c*p*q*R; gfp=(LTFAT_REAL*)gf; for (int r=0;r<c;r++) { for (int w=0;w<R;w++) { for (int l=0;l<q;l++) { for (int k=0;k<p;k++) { negrem = positiverem(k*M-l*a,L); for (s=0;s<d;s++) { rem = (negrem+s*p*M)%L; sbuf[s] = sqrtM*g[r+rem+L*w]; } LTFAT_FFTW(execute)(p_before); for (s=0;s<d2;s++) { gfp[s*ld3] = cbuf[s][0]; gfp[s*ld3+1]= cbuf[s][1]; } gfp+=2; } } } } ltfat_free(sbuf); }