1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206
|
#include "ltfat.h"
#include "ltfat_types.h"
LTFAT_EXTERN void
LTFAT_NAME_COMPLEX(iwfac)(const LTFAT_COMPLEX *gf, const ltfatInt L, const ltfatInt R,
const ltfatInt a, const ltfatInt M, LTFAT_COMPLEX *g)
{
ltfatInt h_a, h_m;
ltfatInt rem, negrem;
LTFAT_REAL scaling, *sbuf, *gfp;
LTFAT_FFTW(plan) p_before;
const ltfatInt b=L/M;
const ltfatInt c=gcd(a, M,&h_a, &h_m);
const ltfatInt p=a/c;
const ltfatInt q=M/c;
const ltfatInt d=b/p;
/* division by d is because of the way FFTW normalizes the transform. */
scaling=1.0/sqrt(M)/d;
sbuf = ltfat_malloc(2*d*sizeof*sbuf);
/* Create plan. In-place. */
p_before = LTFAT_FFTW(plan_dft_1d)(d, (LTFAT_COMPLEX*)sbuf, (LTFAT_COMPLEX*)sbuf,
FFTW_BACKWARD, FFTW_MEASURE);
const ltfatInt ld3=c*p*q*R;
gfp=(LTFAT_REAL*)gf;
for (ltfatInt r=0; r<c; r++)
{
for (ltfatInt w=0; w<R; w++)
{
for (ltfatInt l=0; l<q; l++)
{
for (ltfatInt k=0; k<p; k++)
{
negrem=positiverem(k*M-l*a,L);
for (ltfatInt s=0; s<2*d; s+=2)
{
sbuf[s] = gfp[s*ld3]*scaling;
sbuf[s+1] = gfp[s*ld3+1]*scaling;
}
LTFAT_FFTW(execute)(p_before);
for (ltfatInt s=0; s<d; s++)
{
rem = (negrem+s*p*M)%L;
LTFAT_REAL* gTmp = (LTFAT_REAL*) &(g[r+rem+L*w]);
gTmp[0] = sbuf[2*s];
gTmp[1] = sbuf[2*s+1];
}
gfp+=2;
}
}
}
}
/* Clear the work-array. */
ltfat_free(sbuf);
LTFAT_FFTW(destroy_plan)(p_before);
}
LTFAT_EXTERN void
LTFAT_NAME_REAL(iwfac)(const LTFAT_COMPLEX *gf, const ltfatInt L, const ltfatInt R,
const ltfatInt a, const ltfatInt M, LTFAT_REAL *g)
{
ltfatInt h_a, h_m;
ltfatInt rem, negrem;
LTFAT_REAL scaling, *sbuf, *gfp;
LTFAT_FFTW(plan) p_before;
const ltfatInt b=L/M;
const ltfatInt c=gcd(a, M,&h_a, &h_m);
const ltfatInt p=a/c;
const ltfatInt q=M/c;
const ltfatInt d=b/p;
/* division by d is because of the way FFTW normalizes the transform. */
scaling=1.0/sqrt(M)/d;
sbuf = ltfat_malloc(2*d*sizeof*sbuf);
/* Create plan. In-place. */
p_before = LTFAT_FFTW(plan_dft_1d)(d, (LTFAT_COMPLEX*)sbuf, (LTFAT_COMPLEX*)sbuf,
FFTW_BACKWARD, FFTW_MEASURE);
const ltfatInt ld3=c*p*q*R;
gfp=(LTFAT_REAL*)gf;
for (ltfatInt r=0; r<c; r++)
{
for (ltfatInt w=0; w<R; w++)
{
for (ltfatInt l=0; l<q; l++)
{
for (ltfatInt k=0; k<p; k++)
{
negrem=positiverem(k*M-l*a,L);
for (ltfatInt s=0; s<2*d; s+=2)
{
sbuf[s] = gfp[s*ld3]*scaling;
sbuf[s+1] = gfp[s*ld3+1]*scaling;
}
LTFAT_FFTW(execute)(p_before);
for (ltfatInt s=0; s<d; s++)
{
rem = (negrem+s*p*M)%L;
g[r+rem+L*w] = sbuf[2*s];
}
gfp+=2;
}
}
}
}
/* Clear the work-array. */
ltfat_free(sbuf);
LTFAT_FFTW(destroy_plan)(p_before);
}
LTFAT_EXTERN void
LTFAT_NAME(iwfacreal)(const LTFAT_COMPLEX *gf, const ltfatInt L, const ltfatInt R,
const ltfatInt a, const ltfatInt M, LTFAT_REAL *g)
{
ltfatInt h_a, h_m;
LTFAT_FFTW(plan) p_before;
const ltfatInt b=L/M;
const ltfatInt c=gcd(a, M,&h_a, &h_m);
const ltfatInt p=a/c;
const ltfatInt q=M/c;
const ltfatInt d=b/p;
/* This is a floor operation. */
const ltfatInt d2= d/2+1;
/* division by d is because of the way FFTW normalizes the transform. */
const LTFAT_REAL scaling=1.0/sqrt(M)/d;
LTFAT_REAL *sbuf = ltfat_malloc( d*sizeof*sbuf);
LTFAT_COMPLEX *cbuf = ltfat_malloc(d2*sizeof*cbuf);
/* Create plan. In-place. */
p_before = LTFAT_FFTW(plan_dft_c2r_1d)(d, cbuf, sbuf, FFTW_MEASURE);
const ltfatInt ld3=c*p*q*R;
/* Advancing pointer: Runs through array pointing out the base for the strided operations. */
const LTFAT_COMPLEX *gfp = gf;
for (ltfatInt r=0; r<c; r++)
{
for (ltfatInt w=0; w<R; w++)
{
for (ltfatInt l=0; l<q; l++)
{
for (ltfatInt k=0; k<p; k++)
{
const ltfatInt negrem=positiverem(k*M-l*a,L);
for (ltfatInt s=0; s<d2; s++)
{
cbuf[s] = gfp[s*ld3]*scaling;
}
LTFAT_FFTW(execute)(p_before);
for (ltfatInt s=0; s<d; s++)
{
g[r+(negrem+s*p*M)%L+L*w] = sbuf[s];
}
gfp++;
}
}
}
}
/* Clear the work-arrays. */
LTFAT_SAFEFREEALL(cbuf,sbuf);
LTFAT_FFTW(destroy_plan)(p_before);
}
|