void Forward_SO3_Naive_fftw( int bw, fftw_complex *data, fftw_complex *coeffs, fftw_complex *workspace_cx, fftw_complex *workspace_cx2, double *workspace_re, double *weights, fftw_plan *p1, int flag ) { int j, n, n3 ; int m1, m2 ; int sampHere, coefHere ; int coefHere2 ; int tmpInt ; double *sinPts, *cosPts, *sinPts2, *cosPts2 ; double *wigners, *scratch ; double fudge ; fftw_complex *coeffsPtr ; fftw_complex *dataPtr ; double dn ; n = 2 * bw ; n3 = n * n * n ; /* I'll need these for later */ coeffsPtr = coeffs ; dataPtr = data ; sinPts = workspace_re ; cosPts = sinPts + n ; sinPts2 = cosPts + n ; cosPts2 = sinPts2 + n ; wigners = cosPts2 + n ; scratch = wigners + ( bw * n ) ; /* wigners need at most bw*n space AT ANY given orders m1, m2 */ /* before going further, let's precompute all the sines and cosines I'll need. No matter what order transform I'm doing, these'll stay the same. */ SinEvalPts( n, sinPts ); CosEvalPts( n, cosPts ); SinEvalPts2( n, sinPts2 ); CosEvalPts2( n, cosPts2 ); /* I also need to copy the contents of data to workspace_cx2, given that's where the fftw plan expects data to be. I wish I didn't have to waste so much memory */ memcpy( workspace_cx2, data, sizeof(fftw_complex) * n3 ); /* Stage 1: FFT the "rows". Instead of treating the signal as 3-D object, I can also think of it as an array of size (n^2) x n. This means all I'm doing in the first stage is taking n^2-many FFTs, each of length n. NOTE: Since I'm reusing the FFT code from SpharmonicKit, even though I'm doing the FORWARD SO(3) transform here, I need to call grid_invfourier -> the signs on the complex exponentials are switched (detailed explanation to be put here eventually, but trust me) */ fftw_execute( *p1 ) ; /* normalize the Fourier coefficients (sorry, have to do it) */ /* no, I don't! I can wait till the end */ /* dn = 1. / sqrt( (double) n ) ; for ( j = 0 ; j < n*n*n; j++ ) { workspace_cx[ j ][0] *= dn ; workspace_cx[ j ][1] *= dn ; } */ /* Stage 2: transpose! */ transpose_cx( workspace_cx, workspace_cx2, n*n, n ) ; /* Stage 3: FFT again. */ fftw_execute( *p1 ) ; /* normalize the Fourier coefficients (sorry, have to do it) */ /* no, I don't! I can wait till the end */ /* dn = 1. / ((double) n) ; for ( j = 0 ; j < n*n*n; j++ ) { workspace_cx[ j ][0] *= dn ; workspace_cx[ j ][1] *= dn ; } */ /* Stage 4: transpose again! And note I'm using the tmp space of t2r, t2i again. */ transpose_cx( workspace_cx, workspace_cx2, n*n, n ) ; /* Stage 5: Do the Wigner transforms. This is the tricky bit. Since I'm working with two order indeces, m1 and m2, the for-loops will be more intricate than in the case of the "ordinary" spherical transform on S^2. Also, I will be taking advantage of the symmetries of the Wigner-d functions. As long as I keep my signs and flips right, the Wigner-d's I precompute for an order (m1, m2) transform can generally be used in seven more transforms: (m1,-m2), (m2,m1), (m2,-m1), (-m2,m1), (-m2,-m1), (-m1,m2) and (-m1,-m2). I say "general" because, of course, I'll be transforming at orders (m1,m1), (m1,0) and (0,m1), so I won't get such a huge reduction. Still, this should save time. If assumptions are made regarding the original input signal, e.g. it's strictly real, then one may take advantage of symmetries of the big D wigners (i.e. function of all 3 parameters alpha, beta, gamma) and so simplify the for-loops some and hence increase the speed of the program. However, the for-loops to follow will make no such assumptions. Whether the signal is real or complex, these for-loops will handle it. The for-loops will be "designed" as follows. They will be divided into cases according to the orders: 0) {f_{0,0}} 1) for 0 <= m1 <= bw-1 compute the coefficients i) {f_{ m1, m1}} ii) {f_{-m1,-m1}} iii) {f_{-m1, m1}} iv) {f_{ m1,-m1}} 2) for 1 <= m1 <= bw-1 compute the coefficients i) {f_{ m1, 0}} ii) {f_{-m1, 0}} iii) {f_{ 0, m1}} iv) {f_{ 0,-m1}} 3) for 1 <= m1 <= bw-1 for m1+1 <= m2 <= bw-1 compute the coefficients i) {f_{ m1, m2}} ii) {f_{-m1,-m2}} iii) {f_{ m1,-m2}} iv) {f_{-m1, m2}} v) {f_{ m2, m1}} vi) {f_{-m2,-m1}} vii) {f_{ m2,-m1}} viii) {f_{-m2, m1}} Fasten your seatbelt, folks. It's going to be a bumpy ride. */ /***************************/ /* */ /* {f_{0,0}} coefficient */ /* */ /***************************/ /* compute the wigners I'll need */ genWig_L2( 0, 0, bw, sinPts, cosPts, sinPts2, cosPts2, wigners, scratch ) ; /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( 0, 0, bw ) ; coefHere = coefLoc_so3( 0, 0, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( 0, 0, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; /*** 0 <= m1 <= bw-1 ***/ for ( m1 = 1 ; m1 < bw ; m1 ++ ) { /* compute the wigners I'll need */ genWig_L2( m1, m1, bw, sinPts, cosPts, sinPts2, cosPts2, wigners, scratch ) ; /***************************/ /* */ /* {f_{m1,m1}} coefficient */ /* */ /***************************/ /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, m1, bw ) ; coefHere = coefLoc_so3( m1, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( m1, m1, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; /*****************************/ /* */ /* {f_{-m1,-m1}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* if data is complex */ { /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, -m1, bw ) ; coefHere = coefLoc_so3( -m1, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( -m1, -m1, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m1, m1, bw ) ; coefHere2 = coefLoc_so3( -m1, -m1, bw ) ; for ( j = 0 ; j < bw - m1 ; j ++ ) { coeffs[coefHere2+j][0] = coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -coeffs[coefHere+j][1]; } } /*****************************/ /* */ /* {f_{-m1,m1}} coefficient */ /* */ /*****************************/ /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, m1, bw ) ; coefHere = coefLoc_so3( -m1, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( -m1, m1, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; /*****************************/ /* */ /* {f_{m1,-m1}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* data is complex */ { /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, -m1, bw ) ; coefHere = coefLoc_so3( m1, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( m1, -m1, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( -m1, m1, bw ); coefHere2 = coefLoc_so3( m1, -m1, bw ); for ( j = 0 ; j < bw - m1 ; j ++ ) { coeffs[coefHere2+j][0] = coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -coeffs[coefHere+j][1]; } } } /*** for 1 <= m1 <= bw-1 ***/ for ( m1 = 1 ; m1 < bw ; m1 ++ ) { /* compute the wigners I'll need */ genWig_L2( m1, 0, bw, sinPts, cosPts, sinPts2, cosPts2, wigners, scratch ) ; /***************************/ /* */ /* {f_{m1,0}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, 0, bw ) ; coefHere = coefLoc_so3( m1, 0, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( m1, 0, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; /***************************/ /* */ /* {f_{-m1,0}} coefficient */ /* */ /***************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, 0, bw ) ; coefHere = coefLoc_so3( -m1, 0, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwX( -m1, 0, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m1, 0, bw ); coefHere2 = coefLoc_so3( -m1, 0, bw ); if ( (m1 % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m1 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } /***************************/ /* */ /* {f_{0,m1}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( 0, m1, bw ) ; coefHere = coefLoc_so3( 0, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwX( 0, m1, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; /***************************/ /* */ /* {f_{0,-m1}} coefficient */ /* */ /***************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( 0, -m1, bw ) ; coefHere = coefLoc_so3( 0, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( 0, -m1, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( 0, m1, bw ); coefHere2 = coefLoc_so3( 0, -m1, bw ); if ( (m1 % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m1 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } } /*** 1 <= m1 <= bw-1 m1+1 <= m2 <= bw-1 ***/ for ( m1 = 1 ; m1 < bw ; m1 ++ ) { for ( m2 = m1 + 1 ; m2 < bw ; m2 ++ ) { /* compute the wigners I'll need */ genWig_L2( m1, m2, bw, sinPts, cosPts, sinPts2, cosPts2, wigners, scratch ) ; /***************************/ /* */ /* {f_{m1,m2}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, m2, bw ) ; coefHere = coefLoc_so3( m1, m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( m1, m2, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; /*****************************/ /* */ /* {f_{-m1,-m2}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, -m2, bw ) ; coefHere = coefLoc_so3( -m1, -m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwX( -m1, -m2, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m1, m2, bw ); coefHere2 = coefLoc_so3( -m1, -m2, bw ); if ( ((m2-m1) % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m2 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } /****************************/ /* */ /* {f_{m1,-m2}} coefficient */ /* */ /****************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, -m2, bw ) ; coefHere = coefLoc_so3( m1, -m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( m1, -m2, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; /*****************************/ /* */ /* {f_{-m1,m2}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, m2, bw ) ; coefHere = coefLoc_so3( -m1, m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( -m1, m2, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m1, -m2, bw ); coefHere2 = coefLoc_so3( -m1, m2, bw ); if ( ((m2-m1) % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m2 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } /***************************/ /* */ /* {f_{m2,m1}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m2, m1, bw ) ; coefHere = coefLoc_so3( m2, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwX( m2, m1, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; /*****************************/ /* */ /* {f_{-m2,-m1}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m2, -m1, bw ) ; coefHere = coefLoc_so3( -m2, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( -m2, -m1, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m2, m1, bw ); coefHere2 = coefLoc_so3( -m2, -m1, bw ); if ( ((m2-m1) % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m2 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } /****************************/ /* */ /* {f_{m2,-m1}} coefficient */ /* */ /****************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m2, -m1, bw ) ; coefHere = coefLoc_so3( m2, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( m1, -m2, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; /****************************/ /* */ /* {f_{-m2,m1}} coefficient */ /* */ /****************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m2, m1, bw ) ; coefHere = coefLoc_so3( -m2, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = workspace_cx2 ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( -m1, m2, bw, dataPtr, wigners, weights, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m2, -m1, bw ); coefHere2 = coefLoc_so3( -m2, m1, bw ); if ( ((m2-m1) % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m2 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } } } /* reset coef ptrs */ coeffsPtr = coeffs ; /* need to normalize, one last time */ dn = (M_PI / ( (double) (bw * n )) ); tmpInt = totalCoeffs_so3( bw ) ; for ( j = 0 ; j < tmpInt ; j ++ ) { coeffsPtr[ j ][0] *= dn ; coeffsPtr[ j ][1] *= dn ; } /*** and we're done ! ***/ }
void Inverse_SO3_Naive_fftw( int bw, fftw_complex *coeffs, fftw_complex *data, fftw_complex *workspace_cx, fftw_complex *workspace_cx2, double *workspace_re, fftw_plan *p1, int flag ) { int j, n ; int m1, m2 ; int sampHere , coefHere ; int sampHere2 ; fftw_complex *coeffsPtr, *dataPtr ; double *sinPts, *cosPts, *sinPts2, *cosPts2 ; double *wignersTrans, *scratch ; double dn ; n = 2 * bw ; /* I'll need these for later */ dataPtr = workspace_cx ; ; coeffsPtr = coeffs ; sinPts = workspace_re ; cosPts = sinPts + n ; sinPts2 = cosPts + n ; cosPts2 = sinPts2 + n ; wignersTrans = cosPts2 + n ; scratch = wignersTrans + ( bw * n ) ; /* wignersTrans need at most bw*n space AT ANY given orders m1, m2 */ /* before going further, let's precompute all the sines and cosines I'll need. No matter what order transform I'm doing, these'll stay the same. */ SinEvalPts( n, sinPts ); CosEvalPts( n, cosPts ); SinEvalPts2( n, sinPts2 ); CosEvalPts2( n, cosPts2 ); /* Stage 0.5: Need to normalize the numbers before doing the IDWT */ /* no! I can wait till the end */ /* dn = ( ((double) bw) / M_PI ) ; for ( j = 0 ; j < totalCoeffs_so3( bw ) ; j++ ) { workspace_cx[ j ][0] = coeffs[j][0] * dn ; workspace_cx[ j ][1] = coeffs[j][1] * dn ; } */ /* Stage 1: Do the Inverse Wigner transform. The rcoeffs, icoeffs arrays are assumed to be in the same "arrangement" as that produced by Forward_SO3_Naive(). Since I'm working with two order indeces, m1 and m2, the for-loops will be more intricate than in the case of the "ordinary" spherical transform on S^2. Also, I will be taking advantage of the symmetries of the Wigner-d functions. As long as I keep my signs and flips right, the Wigner-d's I precompute for an order (m1, m2) transform can generally be used in seven more transforms: (m1,-m2), (m2,m1), (m2,-m1), (-m2,m1), (-m2,-m1), (-m1,m2) and (-m1,-m2). The for-loops will be "designed" as follows. They will be divided into cases according to the orders: 0) {f_{0,0}} inverse transform 1) for 0 <= m1 <= bw-1 compute inverse transform of i) {f_{ m1, m1}} ii) {f_{-m1,-m1}} iii) {f_{-m1, m1}} iv) {f_{ m1,-m1}} 2) for 1 <= m1 <= bw-1 compute inverse transform of i) {f_{ m1, 0}} ii) {f_{-m1, 0}} iii) {f_{ 0, m1}} iv) {f_{ 0,-m1}} 3) for 1 <= m1 <= bw-1 for m1+1 <= m2 <= bw-1 compute inverse transform i) {f_{ m1, m2}} ii) {f_{-m1,-m2}} iii) {f_{ m1,-m2}} iv) {f_{-m1, m2}} v) {f_{ m2, m1}} vi) {f_{-m2,-m1}} vii) {f_{ m2,-m1}} viii) {f_{-m2, m1}} If assumptions are made regarding the original input signal, e.g. it's strictly real, then one may take advantage of symmetries of the big D wigners (i.e. function of all 3 parameters alpha, beta, gamma) and so simplify the for-loops some and hence increase the speed of the program. However, the for-loops to follow will make no such assumptions. Whether the signal is real or complex, these for-loops will handle it. Fasten your seatbelt, folks. It's going to be a bumpy ride. */ /* NOTE that I'm using the rdata, idata arrays as tmp space in the early going of the function */ /***************************/ /* */ /* {f_{0,0}} coefficient */ /* */ /***************************/ /* compute the wigners I'll need */ genWigTrans_L2( 0, 0, bw, sinPts, cosPts, sinPts2, cosPts2, wignersTrans, scratch ) ; /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( 0, 0, bw ) ; coefHere = coefLoc_so3( 0, 0, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftw( 0, 0, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; /*** 0 <= m1 <= bw-1 ***/ for ( m1 = 1 ; m1 < bw ; m1 ++ ) { /* compute the wigners I'll need */ genWigTrans_L2( m1, m1, bw, sinPts, cosPts, sinPts2, cosPts2, wignersTrans, scratch ) ; /***************************/ /* */ /* {f_{m1,m1}} coefficient */ /* */ /***************************/ /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, m1, bw ) ; coefHere = coefLoc_so3( m1, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftw( m1, m1, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; /*****************************/ /* */ /* {f_{-m1,-m1}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* if data is complex */ { /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, -m1, bw ) ; coefHere = coefLoc_so3( -m1, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftw( -m1, -m1, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; } else /* otherwise, use symmetry */ { sampHere = sampLoc_so3( m1, m1, bw ); sampHere2 = sampLoc_so3( -m1, -m1, bw ); for ( j = 0 ; j < 2*bw ; j ++ ) { data[sampHere2+j][0] = data[sampHere+j][0]; data[sampHere2+j][1] = -data[sampHere+j][1]; } } /*****************************/ /* */ /* {f_{-m1,m1}} coefficient */ /* */ /*****************************/ /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, m1, bw ) ; coefHere = coefLoc_so3( -m1, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftwY( -m1, m1, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; /*****************************/ /* */ /* {f_{m1,-m1}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* if data is complex */ { /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, -m1, bw ) ; coefHere = coefLoc_so3( m1, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftwY( m1, -m1, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; } else /* otherwise, use symmetry */ { sampHere = sampLoc_so3( -m1, m1, bw ); sampHere2 = sampLoc_so3( m1, -m1, bw ); for ( j = 0 ; j < 2*bw ; j ++ ) { data[sampHere2+j][0] = data[sampHere+j][0]; data[sampHere2+j][1] = -data[sampHere+j][1]; } } } /*** for 1 <= m1 <= bw-1 ***/ for ( m1 = 1 ; m1 < bw ; m1 ++ ) { /* compute the wigners I'll need */ genWigTrans_L2( m1, 0, bw, sinPts, cosPts, sinPts2, cosPts2, wignersTrans, scratch ) ; /***************************/ /* */ /* {f_{m1,0}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, 0, bw ) ; coefHere = coefLoc_so3( m1, 0, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftw( m1, 0, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; /***************************/ /* */ /* {f_{-m1,0}} coefficient */ /* */ /***************************/ if ( flag == 0 ) /* if data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, 0, bw ) ; coefHere = coefLoc_so3( -m1, 0, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftwX( -m1, 0, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; } else /* otherwise, use symmetry */ { sampHere = sampLoc_so3( m1, 0, bw ); sampHere2 = sampLoc_so3( -m1, 0, bw ); for ( j = 0 ; j < 2*bw ; j ++ ) { data[sampHere2+j][0] = data[sampHere+j][0]; data[sampHere2+j][1] = -data[sampHere+j][1]; } } /***************************/ /* */ /* {f_{0,m1}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( 0, m1, bw ) ; coefHere = coefLoc_so3( 0, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftwX( 0, m1, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; /***************************/ /* */ /* {f_{0,-m1}} coefficient */ /* */ /***************************/ if ( flag == 0 ) /* if data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( 0, -m1, bw ) ; coefHere = coefLoc_so3( 0, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftw( 0, -m1, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; } else /* otherwise, use symmetry */ { sampHere = sampLoc_so3( 0, m1, bw ); sampHere2 = sampLoc_so3( 0, -m1, bw ); for ( j = 0 ; j < 2*bw ; j ++ ) { data[sampHere2+j][0] = data[sampHere+j][0]; data[sampHere2+j][1] = -data[sampHere+j][1]; } } } /*** 1 <= m1 <= bw-1 m1+1 <= m2 <= bw-1 ***/ for ( m1 = 1 ; m1 < bw ; m1 ++ ) { for ( m2 = m1 + 1 ; m2 < bw ; m2 ++ ) { /* compute the wigners I'll need */ genWigTrans_L2( m1, m2, bw, sinPts, cosPts, sinPts2, cosPts2, wignersTrans, scratch ) ; /***************************/ /* */ /* {f_{m1,m2}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, m2, bw ) ; coefHere = coefLoc_so3( m1, m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftw( m1, m2, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; /*****************************/ /* */ /* {f_{-m1,-m2}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* if data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, -m2, bw ) ; coefHere = coefLoc_so3( -m1, -m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftwX( -m1, -m2, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; } else /* otherwise, use symmetry */ { sampHere = sampLoc_so3( m1, m2, bw ); sampHere2 = sampLoc_so3( -m1, -m2, bw ); for ( j = 0 ; j < 2*bw ; j ++ ) { data[sampHere2+j][0] = data[sampHere+j][0]; data[sampHere2+j][1] = -data[sampHere+j][1]; } } /****************************/ /* */ /* {f_{m1,-m2}} coefficient */ /* */ /****************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, -m2, bw ) ; coefHere = coefLoc_so3( m1, -m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftwY( m1, -m2, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; /*****************************/ /* */ /* {f_{-m1,m2}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* if data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, m2, bw ) ; coefHere = coefLoc_so3( -m1, m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftwY( -m1, m2, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; } else /* otherwise, use symmetry */ { sampHere = sampLoc_so3( m1, -m2, bw ); sampHere2 = sampLoc_so3( -m1, m2, bw ); for ( j = 0 ; j < 2*bw ; j ++ ) { data[sampHere2+j][0] = data[sampHere+j][0]; data[sampHere2+j][1] = -data[sampHere+j][1]; } } /***************************/ /* */ /* {f_{m2,m1}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m2, m1, bw ) ; coefHere = coefLoc_so3( m2, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftwX( m2, m1, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; /*****************************/ /* */ /* {f_{-m2,-m1}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* if data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m2, -m1, bw ) ; coefHere = coefLoc_so3( -m2, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftw( -m2, -m1, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; } else /* otherwise, use symmetry */ { sampHere = sampLoc_so3( m2, m1, bw ); sampHere2 = sampLoc_so3( -m2, -m1, bw ); for ( j = 0 ; j < 2*bw ; j ++ ) { data[sampHere2+j][0] = data[sampHere+j][0]; data[sampHere2+j][1] = -data[sampHere+j][1]; } } /****************************/ /* */ /* {f_{m2,-m1}} coefficient */ /* */ /****************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m2, -m1, bw ) ; coefHere = coefLoc_so3( m2, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftwY( m1, -m2, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; /****************************/ /* */ /* {f_{-m2,m1}} coefficient */ /* */ /****************************/ if ( flag == 0 ) /* if data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m2, m1, bw ) ; coefHere = coefLoc_so3( -m2, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveSynthesis_fftwY( -m1, m2, bw, coeffsPtr, wignersTrans, dataPtr, workspace_cx2 ) ; } else /* otherwise, use symmetry */ { sampHere = sampLoc_so3( m2, -m1, bw ); sampHere2 = sampLoc_so3( -m2, m1, bw ); for ( j = 0 ; j < 2*bw ; j ++ ) { data[sampHere2+j][0] = data[sampHere+j][0]; data[sampHere2+j][1] = -data[sampHere+j][1]; } } } } /* I need to set some zeros so that I can take the fft correctly */ /* reset ptrs to correct starting positions */ dataPtr = data + (n)*(bw) ; for ( m1 = 0 ; m1 < bw ; m1 ++ ) { memset( dataPtr, 0, sizeof(fftw_complex) * n ); dataPtr += (2*n)*(bw) ; } dataPtr = data + bw*n*(n); memset( dataPtr, 0, sizeof(fftw_complex) * n * n ); dataPtr += n * n + n*bw; for ( m1 = 1 ; m1 < bw ; m1 ++ ) { memset( dataPtr, 0, sizeof(fftw_complex) * n ); dataPtr += (2*n)*(bw) ; } /* Stage 2: transpose! Note I'm using the rdata, idata arrays as tmp space */ transpose_cx( data, workspace_cx, n, n*n ) ; /* Stage 3: FFT the "rows". Instead of treating the signal as 3-D object, I can also think of it as an array of size (n^2) x n. This means all I'm doing in the first stage is taking n^2-many FFTs, each of length n. NOTE: Since I'm reusing the FFT code from SpharmonicKit, even though I'm doing the INVERSE SO(3) transform here, I need to call grid_fourier -> the signs on the complex exponentials are switched (detailed explanation to be put here eventually, but trust me) */ fftw_execute( *p1 ) ; /* normalize the Fourier coefficients (sorry, have to do it) */ /* no! I can wait till the end */ /* dn = sqrt( (double) n ); for ( j = 0 ; j < n*n*n; j++ ) { data[ j ][0] *= dn ; data[ j ][1] *= dn ; } dn = 1./( (double) n ); for ( j = 0 ; j < n*n*n; j++ ) { data[ j ][0] *= dn ; data[ j ][1] *= dn ; } */ /* Stage 4: transpose! Note I'm using the rdata, idata arrays as tmp space */ transpose_cx( data, workspace_cx, n, n*n ) ; /* Stage 5: FFT again. Note that THIS TIME, the rdata, idata arrays will hold the final answers I want */ fftw_execute( *p1 ) ; /* normalize the Fourier coefficients (sorry, have to do it) */ dn = 1./( (double) n ); dn *= ( ((double) bw) / M_PI ) ; for ( j = 0 ; j < n*n*n; j++ ) { data[ j ][0] *= dn ; data[ j ][1] *= dn ; } /* and that's all, folks */ }
int main ( int argc , char **argv ) { int i, j, m1, m2, bw, n ; int loops, m ; long int seed ; double *coeffs, *signal, *newcoeffs; double *wigners, *wignersTrans ; double *workspace, *scratch ; double *weights ; double *sinPts, *cosPts ; double *sinPts2, *cosPts2 ; double tmp_error, sum_error; double tmp_relerror, sum_relerror; double tstartA, tstopA, runtimeA ; double tstartB, tstopB, runtimeB ; double *relerror, *curmax; double ave_error, ave_relerror, stddev_error, stddev_relerror ; FILE *fp ; if (argc < 5) { fprintf(stdout,"Usage: test_Wigner_Naive m1 m2 bw loops [output_file]\n"); exit(0); } m1 = atoi( argv[1] ); m2 = atoi( argv[2] ); bw = atoi( argv[3] ); loops = atoi( argv[4] ) ; m = MAX( ABS( m1 ) , ABS( m2 ) ) ; n = 2 * bw ; runtimeA = 0.0 ; runtimeB = 0.0 ; weights = ( double * ) malloc(sizeof( double ) * (2*bw) ) ; coeffs = ( double * ) malloc(sizeof( double ) * (bw - m) ) ; newcoeffs = ( double * ) malloc(sizeof( double ) * (bw - m) ) ; signal = ( double * ) malloc(sizeof( double ) * n ) ; wigners = ( double * ) malloc( sizeof( double ) * ( bw - m ) * n ) ; wignersTrans = ( double * ) malloc( sizeof( double ) * ( bw - m ) * n ) ; workspace = (double *) malloc(sizeof( double ) * (4 + 6) * n ) ; sinPts = workspace ; cosPts = sinPts + n ; sinPts2 = cosPts + n ; cosPts2 = sinPts2 + n ; scratch = cosPts2 + n ; /* scratch needs to be of size 6*n */ /* note that the definition of wigSpec requires that instead of evaluating at beta, I need to evaluate at beta/2; ergo I call SinEvalPts2 instead of SinEvalPts, etc etc */ /* generate seed for random number generator */ time ( &seed ) ; srand48( seed ) ; /* precompute sines and cosines appropriate for making the wigners */ SinEvalPts( n, sinPts ) ; CosEvalPts( n, cosPts ) ; SinEvalPts2( n, sinPts2 ) ; CosEvalPts2( n, cosPts2 ) ; /* make quadrature weights */ makeweights2( bw, weights ); /* make the wigners */ genWig_L2( m1, m2, bw, sinPts, cosPts, sinPts2, cosPts2, wigners, scratch ) ; /* now make the wigners - transpose version! */ genWigTrans_L2( m1, m2, bw, sinPts, cosPts, sinPts2, cosPts2, wignersTrans, scratch ) ; /** space for errors **/ relerror = (double *) malloc(sizeof(double) * loops); curmax = (double *) malloc(sizeof(double) * loops); sum_error = 0.0 ; sum_relerror = 0.0 ; for ( i = 0 ; i < loops ; i ++ ) { /* generate random coeffs */ for( j = 0 ; j < (bw - m) ; j++ ) coeffs[ j ] = drand48() ; /* turn on stop watch */ tstartA = csecond () ; /* now synthesize */ wigNaiveSynthesis( m1, m2, bw, coeffs, wignersTrans, signal, scratch ) ; tstopA = csecond () ; runtimeA += (tstopA - tstartA); tstartB = csecond () ; /* now analyze */ wigNaiveAnalysis( m1, m2, bw, signal, wigners, weights, newcoeffs, scratch ) ; /* turn off stop watch */ tstopB = csecond () ; runtimeB += (tstopB - tstartB); relerror[ i ] = 0.0 ; curmax[ i ] = 0.0 ; /* now figure out errors */ for( j = 0 ; j < bw - m ; j ++ ) { tmp_error = fabs( coeffs[j] - newcoeffs[j] ); tmp_relerror = tmp_error / ( fabs( coeffs[j] ) + pow( 10.0, -50.0 ) ); curmax[ i ] = MAX( curmax[ i ], tmp_error ); relerror[ i ] = MAX( relerror[ i ], tmp_relerror ); } sum_error += curmax[ i ] ; sum_relerror += relerror[ i ] ; } ave_error = sum_error / ( (double) loops ); ave_relerror = sum_relerror / ( (double) loops ); stddev_error = 0.0 ; stddev_relerror = 0.0; for( i = 0 ; i < loops ; i ++ ) { stddev_error += pow( ave_error - curmax[ i ] , 2.0 ); stddev_relerror += pow( ave_relerror - relerror[ i ] , 2.0 ); } /*** this won't work if loops == 1 ***/ if( loops != 1 ) { stddev_error = sqrt(stddev_error / ( (double) (loops - 1) ) ); stddev_relerror = sqrt(stddev_relerror / ( (double) (loops - 1) ) ); } fprintf(stderr,"bw = %d\tm1 = %d\tm2 = %d\n",bw, m1, m2); fprintf(stderr,"total runtime: %.4e seconds\n", runtimeA+runtimeB); fprintf(stderr,"average forward runtime: %.4e seconds per iteration\n", runtimeB/((double) loops)); fprintf(stderr,"average inverse runtime: %.4e seconds per iteration\n", runtimeA/((double) loops)); fprintf(stderr,"Average r-o error:\t\t %.4e\t", sum_error/((double) loops)); fprintf(stderr,"std dev: %.4e\n",stddev_error); fprintf(stderr,"Average (r-o)/o error:\t\t %.4e\t", sum_relerror/((double) loops)); fprintf(stderr,"std dev: %.4e\n\n",stddev_relerror); if ( argc == 6 ) { fp = fopen(argv[5], "w"); for ( i = 0 ; i < bw - m ; i ++ ) fprintf(fp,"%.16f\n", coeffs[i] - newcoeffs[i]); } free( curmax ) ; free( relerror ) ; free( workspace ) ; free( wignersTrans ) ; free( wigners ) ; free( signal ) ; free( newcoeffs ) ; free( coeffs ) ; free( weights ) ; return 0 ; }
void Forward_SO3_Naive_fftw_wo( int bw, fftw_complex *data, fftw_complex *coeffs, fftw_complex *workspace_cx, REAL *workspace_re, fftw_plan *p1, int flag ) { int j, n ; int m1, m2 ; int sampHere, coefHere ; int coefHere2 ; double *sinPts, *cosPts, *sinPts2, *cosPts2 ; double *wigners, *scratch ; double fudge ; fftw_complex *coeffsPtr ; fftw_complex *dataPtr ; double dn ; n = 2 * bw ; sinPts = workspace_re ; cosPts = sinPts + n ; sinPts2 = cosPts + n ; cosPts2 = sinPts2 + n ; wigners = cosPts2 + n ; scratch = wigners + ( bw * n ) ; /* wigners need at most bw*n space AT ANY given orders m1, m2 */ /* before going further, let's precompute all the sines and cosines I'll need. No matter what order transform I'm doing, these'll stay the same. */ SinEvalPts( n, sinPts ); CosEvalPts( n, cosPts ); SinEvalPts2( n, sinPts2 ); CosEvalPts2( n, cosPts2 ); /* Stage 1: FFT */ fftw_execute( *p1 ) ; /* Stage 2: Do the Wigner transforms. This is the tricky bit. Since I'm working with two order indeces, m1 and m2, the for-loops will be more intricate than in the case of the "ordinary" spherical transform on S^2. Also, I will be taking advantage of the symmetries of the Wigner-d functions. As long as I keep my signs and flips right, the Wigner-d's I precompute for an order (m1, m2) transform can generally be used in seven more transforms: (m1,-m2), (m2,m1), (m2,-m1), (-m2,m1), (-m2,-m1), (-m1,m2) and (-m1,-m2). I say "general" because, of course, I'll be transforming at orders (m1,m1), (m1,0) and (0,m1), so I won't get such a huge reduction. Still, this should save time. If assumptions are made regarding the original input signal, e.g. it's strictly real, then one may take advantage of symmetries of the big D wigners (i.e. function of all 3 parameters alpha, beta, gamma) and so simplify the for-loops some and hence increase the speed of the program. However, the for-loops to follow will make no such assumptions. Whether the signal is real or complex, these for-loops will handle it. The for-loops will be "designed" as follows. They will be divided into cases according to the orders: 0) {f_{0,0}} 1) for 0 <= m1 <= bw-1 compute the coefficients i) {f_{ m1, m1}} ii) {f_{-m1,-m1}} iii) {f_{-m1, m1}} iv) {f_{ m1,-m1}} 2) for 1 <= m1 <= bw-1 compute the coefficients i) {f_{ m1, 0}} ii) {f_{-m1, 0}} iii) {f_{ 0, m1}} iv) {f_{ 0,-m1}} 3) for 1 <= m1 <= bw-1 for m1+1 <= m2 <= bw-1 compute the coefficients i) {f_{ m1, m2}} ii) {f_{-m1,-m2}} iii) {f_{ m1,-m2}} iv) {f_{-m1, m2}} v) {f_{ m2, m1}} vi) {f_{-m2,-m1}} vii) {f_{ m2,-m1}} viii) {f_{-m2, m1}} Fasten your seatbelt, folks. It's going to be a bumpy ride. */ /***************************/ /* */ /* {f_{0,0}} coefficient */ /* */ /***************************/ /* compute the wigners I'll need */ genWig_L2( 0, 0, bw, sinPts, cosPts, sinPts2, cosPts2, wigners, scratch ) ; /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( 0, 0, bw ) ; coefHere = coefLoc_so3( 0, 0, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( 0, 0, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; /*** 0 <= m1 <= bw-1 ***/ for ( m1 = 1 ; m1 < bw ; m1 ++ ) { /* compute the wigners I'll need */ genWig_L2( m1, m1, bw, sinPts, cosPts, sinPts2, cosPts2, wigners, scratch ) ; /***************************/ /* */ /* {f_{m1,m1}} coefficient */ /* */ /***************************/ /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, m1, bw ) ; coefHere = coefLoc_so3( m1, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( m1, m1, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; /*****************************/ /* */ /* {f_{-m1,-m1}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* if data is complex */ { /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, -m1, bw ) ; coefHere = coefLoc_so3( -m1, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( -m1, -m1, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m1, m1, bw ) ; coefHere2 = coefLoc_so3( -m1, -m1, bw ) ; for ( j = 0 ; j < bw - m1 ; j ++ ) { coeffs[coefHere2+j][0] = coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -coeffs[coefHere+j][1]; } } /*****************************/ /* */ /* {f_{-m1,m1}} coefficient */ /* */ /*****************************/ /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, m1, bw ) ; coefHere = coefLoc_so3( -m1, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( -m1, m1, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; /*****************************/ /* */ /* {f_{m1,-m1}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* data is complex */ { /* now, get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, -m1, bw ) ; coefHere = coefLoc_so3( m1, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( m1, -m1, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( -m1, m1, bw ); coefHere2 = coefLoc_so3( m1, -m1, bw ); for ( j = 0 ; j < bw - m1 ; j ++ ) { coeffs[coefHere2+j][0] = coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -1.*coeffs[coefHere+j][1]; } } } /*** for 1 <= m1 <= bw-1 ***/ for ( m1 = 1 ; m1 < bw ; m1 ++ ) { /* compute the wigners I'll need */ genWig_L2( m1, 0, bw, sinPts, cosPts, sinPts2, cosPts2, wigners, scratch ) ; /***************************/ /* */ /* {f_{m1,0}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, 0, bw ) ; coefHere = coefLoc_so3( m1, 0, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( m1, 0, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; /***************************/ /* */ /* {f_{-m1,0}} coefficient */ /* */ /***************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, 0, bw ) ; coefHere = coefLoc_so3( -m1, 0, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwX( -m1, 0, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m1, 0, bw ); coefHere2 = coefLoc_so3( -m1, 0, bw ); if ( (m1 % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m1 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } /***************************/ /* */ /* {f_{0,m1}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( 0, m1, bw ) ; coefHere = coefLoc_so3( 0, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwX( 0, m1, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; /***************************/ /* */ /* {f_{0,-m1}} coefficient */ /* */ /***************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( 0, -m1, bw ) ; coefHere = coefLoc_so3( 0, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( 0, -m1, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( 0, m1, bw ); coefHere2 = coefLoc_so3( 0, -m1, bw ); if ( (m1 % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m1 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } } /*** 1 <= m1 <= bw-1 m1+1 <= m2 <= bw-1 ***/ for ( m1 = 1 ; m1 < bw ; m1 ++ ) { for ( m2 = m1 + 1 ; m2 < bw ; m2 ++ ) { /* compute the wigners I'll need */ genWig_L2( m1, m2, bw, sinPts, cosPts, sinPts2, cosPts2, wigners, scratch ) ; /***************************/ /* */ /* {f_{m1,m2}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, m2, bw ) ; coefHere = coefLoc_so3( m1, m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( m1, m2, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; /*****************************/ /* */ /* {f_{-m1,-m2}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, -m2, bw ) ; coefHere = coefLoc_so3( -m1, -m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwX( -m1, -m2, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m1, m2, bw ); coefHere2 = coefLoc_so3( -m1, -m2, bw ); if ( ((m2-m1) % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m2 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } /****************************/ /* */ /* {f_{m1,-m2}} coefficient */ /* */ /****************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m1, -m2, bw ) ; coefHere = coefLoc_so3( m1, -m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( m1, -m2, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; /*****************************/ /* */ /* {f_{-m1,m2}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m1, m2, bw ) ; coefHere = coefLoc_so3( -m1, m2, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( -m1, m2, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m1, -m2, bw ); coefHere2 = coefLoc_so3( -m1, m2, bw ); if ( ((m2-m1) % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m2 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } /***************************/ /* */ /* {f_{m2,m1}} coefficient */ /* */ /***************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m2, m1, bw ) ; coefHere = coefLoc_so3( m2, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwX( m2, m1, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; /*****************************/ /* */ /* {f_{-m2,-m1}} coefficient */ /* */ /*****************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m2, -m1, bw ) ; coefHere = coefLoc_so3( -m2, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftw( -m2, -m1, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m2, m1, bw ); coefHere2 = coefLoc_so3( -m2, -m1, bw ); if ( ((m2-m1) % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m2 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } /****************************/ /* */ /* {f_{m2,-m1}} coefficient */ /* */ /****************************/ /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( m2, -m1, bw ) ; coefHere = coefLoc_so3( m2, -m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( m1, -m2, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; /****************************/ /* */ /* {f_{-m2,m1}} coefficient */ /* */ /****************************/ if ( flag == 0 ) /* data is complex */ { /* get the locations of where the samples I have to transform are, and where the coefficients have to go */ sampHere = sampLoc_so3( -m2, m1, bw ) ; coefHere = coefLoc_so3( -m2, m1, bw ) ; /* ok, reset sample, coef ptrs */ coeffsPtr = coeffs ; dataPtr = data ; /* now advance by the computed amounts */ dataPtr += sampHere ; coeffsPtr += coefHere ; /* now transform the real and imaginary parts of the data */ wigNaiveAnalysis_fftwY( -m1, m2, bw, dataPtr, wigners, coeffsPtr, workspace_cx ) ; } else /* data is real, so use symmetry */ { coefHere = coefLoc_so3( m2, -m1, bw ); coefHere2 = coefLoc_so3( -m2, m1, bw ); if ( ((m2-m1) % 2) == 0 ) fudge = 1.0 ; else fudge = -1.0 ; for ( j = 0 ; j < bw - m2 ; j ++ ) { coeffs[coefHere2+j][0] = fudge * coeffs[coefHere+j][0]; coeffs[coefHere2+j][1] = -fudge * coeffs[coefHere+j][1]; } } } } /* reset coef ptrs */ coeffsPtr = coeffs ; /* need to normalize, one last time */ dn = (M_PI / ( (double) (bw * n )) ); for ( j = 0 ; j < totalCoeffs_so3( bw ) ; j ++ ) { coeffsPtr[ j ][0] *= dn ; coeffsPtr[ j ][1] *= dn ; } /*** and we're done ! ***/ }
int main ( int argc , char **argv ) { int i, m1, m2, bw, n ; int m ; double *workspace, *scratch ; double *sinPts, *cosPts, *result ; double *sinPts2, *cosPts2 ; FILE *fp ; if (argc < 5) { fprintf(stdout,"Usage: test_genWig m1 m2 bw output_file_name\n"); exit(0); } m1 = atoi( argv[1] ); m2 = atoi( argv[2] ); bw = atoi( argv[3] ); m = MAX( ABS( m1 ) , ABS( m2 ) ) ; n = 2 * bw ; result = ( double * ) malloc(sizeof( double ) * n * ( bw - m ) ) ; workspace = (double *) malloc(sizeof( double ) * (4 + 6) * n ) ; sinPts = workspace ; cosPts = sinPts + n ; sinPts2 = cosPts + n ; cosPts2 = sinPts2 + n ; scratch = cosPts2 + n ; /* scratch needs to be of size 6*n */ /* Compute appropriate sines and cosines at Chebyshev points (or their slight variants) note that the definition of wigSpec requires that instead of evaluating at beta, I need to evaluate at beta/2; ergo I call SinEvalPts2 instead of SinEvalPts, etc etc */ SinEvalPts( n, sinPts ) ; CosEvalPts( n, cosPts ) ; SinEvalPts2( n, sinPts2 ) ; CosEvalPts2( n, cosPts2 ) ; genWig_L2( m1, m2, bw, sinPts, cosPts, sinPts2, cosPts2, result, scratch ) ; fp = fopen( argv[4], "w" ); for ( i = 0 ; i < n*(bw-m) ; i++ ) fprintf(fp, "%.15f\n", result[i]); fclose( fp ) ; free( workspace ) ; free( result ) ; return 0 ; }